]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliPhysicsSelection.cxx
- AliPhysicsSelection: protected writing of fHistBunchCrossing and fHistTriggerPattern
[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 //   Origin: Jan Fiete Grosse-Oetringhaus, CERN 
96 //           Michele Floris, CERN
97 //-------------------------------------------------------------------------
98
99 #include <Riostream.h>
100 #include <TH1F.h>
101 #include <TH2F.h>
102 #include <TList.h>
103 #include <TIterator.h>
104 #include <TDirectory.h>
105 #include <TObjArray.h>
106
107 #include <AliPhysicsSelection.h>
108
109 #include <AliTriggerAnalysis.h>
110 #include <AliLog.h>
111
112 #include <AliESDEvent.h>
113 #include <AliAnalysisTaskSE.h>
114 #include "AliAnalysisManager.h"
115
116 ClassImp(AliPhysicsSelection)
117
118 AliPhysicsSelection::AliPhysicsSelection() :
119   AliAnalysisCuts("AliPhysicsSelection", "AliPhysicsSelection"),
120   fCurrentRun(-1),
121   fMC(kFALSE),
122   fCollTrigClasses(),
123   fBGTrigClasses(),
124   fTriggerAnalysis(),
125   fBackgroundIdentification(0),
126   fHistBunchCrossing(0),
127   fHistTriggerPattern(0),
128   fSkipTriggerClassSelection(0),
129   fUsingCustomClasses(0),
130   fSkipV0(0),
131   fBIFactorA(1),
132   fBIFactorC(1),
133   fRatioBEEE(2),
134   fComputeBG(0),
135   fUseBXNumbers(1),
136   fUseMuonTriggers(0),
137   fFillingScheme(""),
138   fBin0CallBack(""),
139   fBin0CallBackPointer(0)
140 {
141   // constructor
142   
143   fCollTrigClasses.SetOwner(1);
144   fBGTrigClasses.SetOwner(1);
145   fTriggerAnalysis.SetOwner(1);
146   fHistStatistics[0] = 0;
147   fHistStatistics[1] = 0;
148   
149   AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
150 }
151     
152 AliPhysicsSelection::~AliPhysicsSelection()
153 {
154   // destructor
155
156   fCollTrigClasses.Delete();
157   fBGTrigClasses.Delete();
158   fTriggerAnalysis.Delete();
159
160   if (fHistStatistics[0])
161   {
162     delete fHistStatistics[0];
163     fHistStatistics[0] = 0;
164   }
165   if (fHistStatistics[1])
166   {
167     delete fHistStatistics[1];
168     fHistStatistics[1] = 0;
169   }
170   
171   if (fHistBunchCrossing)
172   {
173     delete fHistBunchCrossing;
174     fHistBunchCrossing = 0;
175   }
176   if (fHistTriggerPattern)
177   {
178     delete fHistTriggerPattern;
179     fHistTriggerPattern = 0;
180   }
181
182 }
183
184 Bool_t AliPhysicsSelection::CheckTriggerClass(const AliESDEvent* aEsd, const char* trigger) const
185 {
186   // checks if the given trigger class(es) are found for the current event
187   // format of trigger: +TRIGGER1 -TRIGGER2
188   //   requires TRIGGER1 and rejects TRIGGER2
189   
190   Bool_t foundBCRequirement = kFALSE;
191   Bool_t foundCorrectBC = kFALSE;
192   
193   TString str(trigger);
194   TObjArray* tokens = str.Tokenize(" ");
195   
196   for (Int_t i=0; i < tokens->GetEntries(); i++)
197   {
198     TString str2(((TObjString*) tokens->At(i))->String());
199     
200     if (str2[0] == '+' || str2[0] == '-')
201     {
202       Bool_t flag = (str2[0] == '+');
203       
204       str2.Remove(0, 1);
205       
206       if (flag && !aEsd->IsTriggerClassFired(str2))
207       {
208         AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", str2.Data()));
209         delete tokens;
210         return kFALSE;
211       }
212       if (!flag && aEsd->IsTriggerClassFired(str2))
213       {
214         AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", str2.Data()));
215         delete tokens;
216         return kFALSE;
217       }
218     }
219     else if (str2[0] == '#')
220     {
221       foundBCRequirement = kTRUE;
222     
223       str2.Remove(0, 1);
224       
225       Int_t bcNumber = str2.Atoi();
226       AliDebug(AliLog::kDebug, Form("Checking for bunch crossing number %d", bcNumber));
227       
228       if (aEsd->GetBunchCrossNumber() == bcNumber)
229       {
230         foundCorrectBC = kTRUE;
231         AliDebug(AliLog::kDebug, Form("Found correct bunch crossing %d", bcNumber));
232       }
233     }
234     else
235       AliFatal(Form("Invalid trigger syntax: %s", trigger));
236   }
237   
238   delete tokens;
239   
240   if (foundBCRequirement && !foundCorrectBC)
241     return kFALSE;
242   
243   return kTRUE;
244 }
245     
246 Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
247 {
248   // checks if the given event is a collision candidate
249   
250   if (fCurrentRun != aEsd->GetRunNumber())
251     if (!Initialize(aEsd->GetRunNumber()))
252       AliFatal(Form("Could not initialize for run %d", aEsd->GetRunNumber()));
253     
254   const AliESDHeader* esdHeader = aEsd->GetHeader();
255   if (!esdHeader)
256   {
257     AliError("ESD Header could not be retrieved");
258     return kFALSE;
259   }
260   
261   // check event type; should be PHYSICS = 7 for data and 0 for MC
262   if (!fMC)
263   {
264     if (esdHeader->GetEventType() != 7)
265       return kFALSE;
266   }
267   else
268   {
269     if (esdHeader->GetEventType() != 0)
270       AliFatal(Form("Invalid event type for MC: %d", esdHeader->GetEventType()));
271   }
272   
273   Bool_t accept = kFALSE;
274     
275   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
276   for (Int_t i=0; i < count; i++)
277   {
278     const char* triggerClass = 0;
279     if (i < fCollTrigClasses.GetEntries())
280       triggerClass = ((TObjString*) fCollTrigClasses.At(i))->String();
281     else
282       triggerClass = ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
283   
284     AliDebug(AliLog::kDebug, Form("Processing trigger class %s", triggerClass));
285   
286     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
287   
288     triggerAnalysis->FillTriggerClasses(aEsd);
289     
290     if (CheckTriggerClass(aEsd, triggerClass))
291     {
292       triggerAnalysis->FillHistograms(aEsd);
293   
294       Bool_t isBin0 = kFALSE;
295       if (fBin0CallBack != "") {
296         AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
297         if (!mgr) {
298           AliError("Cannot get the analysis manager");
299         }  
300         else {
301           isBin0 = ((AliAnalysisTaskSE*)mgr->GetTask(fBin0CallBack.Data()))->IsEventInBinZero();
302         }
303       } else if (fBin0CallBackPointer) {
304           isBin0 = (*fBin0CallBackPointer)(aEsd);
305         
306       }
307       
308
309       
310       // hardware trigger (should only remove events for MC)
311       // replay CINT1B hardware trigger
312       // TODO this has to depend on the actual hardware trigger (and that depends on the run...)
313       Int_t fastORHW = triggerAnalysis->SPDFiredChips(aEsd, 1); // SPD number of chips from trigger bits (!)
314       Bool_t v0A       = fSkipV0 ? 0 :triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0A);
315       Bool_t v0C       = fSkipV0 ? 0 :triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0C);
316       Bool_t v0AHW     = fSkipV0 ? 0 :(triggerAnalysis->V0Trigger(aEsd, AliTriggerAnalysis::kASide, kTRUE) == AliTriggerAnalysis::kV0BB);// should replay hw trigger
317       Bool_t v0CHW     = fSkipV0 ? 0 :(triggerAnalysis->V0Trigger(aEsd, AliTriggerAnalysis::kCSide, kTRUE) == AliTriggerAnalysis::kV0BB);// should replay hw trigger
318       // offline trigger
319       Int_t fastOROffline = triggerAnalysis->SPDFiredChips(aEsd, 0); // SPD number of chips from clusters (!)
320       Bool_t v0ABG = fSkipV0 ? 0 :triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0ABG);
321       Bool_t v0CBG = fSkipV0 ? 0 :triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0CBG);
322       Bool_t v0BG = v0ABG || v0CBG;
323
324       // fmd
325       Bool_t fmdA =  triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kFMDA);
326       Bool_t fmdC =  triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kFMDC);
327       Bool_t fmd  = fmdA || fmdC;
328     
329       // SSD
330       Int_t ssdClusters = triggerAnalysis->SSDClusters(aEsd);
331
332       // Some "macros"
333       Bool_t mb1 = (fastOROffline > 0 || v0A || v0C) && (!v0BG);
334       Bool_t mb1prime = (fastOROffline > 1 || (fastOROffline > 0 && (v0A || v0C)) || (v0A && v0C) ) && (!v0BG);
335
336       // Background rejection
337       Bool_t bgID = ! fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd));
338
339       Int_t ntrig = fastOROffline; // any 2 hits
340       if(v0A)              ntrig += 1;
341       if(v0C)              ntrig += 1; //v0C alone is enough
342       if(fmd)              ntrig += 1;
343       if(ssdClusters>1)    ntrig += 1;
344
345
346       Bool_t hwTrig = fastORHW > 0 || v0AHW || v0CHW;
347
348       // Fill trigger pattern histo
349       Int_t tpatt = 0;
350       if (fastORHW>0) tpatt+=1;
351       if (v0AHW)      tpatt+=2;
352       if (v0CHW)      tpatt+=4;
353       fHistTriggerPattern->Fill( tpatt );
354
355       // fill statistics and return decision
356       const Int_t nHistStat = 2;
357       for(Int_t iHistStat = 0; iHistStat < nHistStat; iHistStat++){
358         if (iHistStat == kStatIdxBin0 && !isBin0) continue; // skip the filling of bin0 stats if the event is not in the bin0
359       
360         fHistStatistics[iHistStat]->Fill(kStatTriggerClass, i);
361
362
363
364         // We fill the rest only if hw trigger is ok
365         if (!hwTrig)
366           {
367             AliDebug(AliLog::kDebug, "Rejecting event because hardware trigger is not fired");
368             continue;
369           } else {       
370           fHistStatistics[iHistStat]->Fill(kStatHWTrig, i);
371         }
372       
373
374         // v0 BG stats
375         if (v0ABG)
376           fHistStatistics[iHistStat]->Fill(kStatV0ABG, i);
377         if (v0CBG)
378           fHistStatistics[iHistStat]->Fill(kStatV0CBG, i);
379
380         // We fill the rest only if mb1 && ! v0BG
381         if (mb1)
382           fHistStatistics[iHistStat]->Fill(kStatMB1, i);
383         else continue;
384
385         if (mb1prime)
386           fHistStatistics[iHistStat]->Fill(kStatMB1Prime, i);
387
388         if (fmd)
389           fHistStatistics[iHistStat]->Fill(kStatFMD, i);
390
391         if(ssdClusters>1)
392           fHistStatistics[iHistStat]->Fill(kStatSSD1, i);
393
394         if(ntrig >= 2 && !v0BG) 
395           fHistStatistics[iHistStat]->Fill(kStatAny2Hits, i);
396
397         if (fastOROffline > 0)
398           fHistStatistics[iHistStat]->Fill(kStatFO1, i);
399         if (fastOROffline > 1)
400           fHistStatistics[iHistStat]->Fill(kStatFO2, i);
401         
402         if (v0A)
403           fHistStatistics[iHistStat]->Fill(kStatV0A, i);
404         if (v0C)
405           fHistStatistics[iHistStat]->Fill(kStatV0C, i);
406         
407         //       if (fastOROffline > 1 && !v0BG)
408         //         fHistStatistics[iHistStat]->Fill(kStatFO2NoBG, i);
409             
410         if (fastOROffline > 0 && (v0A || v0C) && !v0BG)
411           fHistStatistics[iHistStat]->Fill(kStatFO1AndV0, i);
412   
413         if (v0A && v0C && !v0BG && !bgID)
414           fHistStatistics[iHistStat]->Fill(kStatV0, i);
415
416
417         if ( mb1 )
418         
419           {
420             if (!v0BG || fSkipV0)
421               {
422                 if (!v0BG) fHistStatistics[iHistStat]->Fill(kStatOffline, i);
423       
424                 if (fBackgroundIdentification && bgID)
425                   {
426                     AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
427                     fHistStatistics[iHistStat]->Fill(kStatBG, i);
428                   }
429                 else
430                   {
431                     AliDebug(AliLog::kDebug, "Accepted event for histograms");
432             
433                     fHistStatistics[iHistStat]->Fill(kStatAccepted, i);
434                     if(iHistStat == kStatIdxAll) fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i); // Fill only for all (avoid double counting)
435                     if((i < fCollTrigClasses.GetEntries() || fSkipTriggerClassSelection) && (iHistStat==kStatIdxAll))
436                       accept = kTRUE; // only set for "all" (should not really matter)
437                   }
438               }
439             else
440               AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
441           }
442         else
443           AliDebug(AliLog::kDebug, "Rejecting event because trigger condition is not fulfilled");
444       }
445     }
446   }
447  
448   if (accept)
449     AliDebug(AliLog::kDebug, "Accepted event as collision candidate");
450   
451   return accept;
452 }
453
454 Int_t AliPhysicsSelection::GetTriggerScheme(UInt_t runNumber) const
455 {
456   // returns the current trigger scheme (classes that are accepted/rejected)
457   
458   if (fMC)
459     return 0;
460     
461   // TODO dependent on run number
462   
463   switch (runNumber)
464   {
465     // CSMBB triggers
466     case 104044:
467     case 105054:
468     case 105057:
469       return 2;
470   }
471   
472   // default: CINT1 suite
473   return 1;
474 }  
475
476 const char * AliPhysicsSelection::GetFillingScheme(UInt_t runNumber)  {
477
478   if(fMC) return "MC";
479
480   if      (runNumber >= 104065 && runNumber <= 104160) {
481     return "4x4a";
482   } 
483   else if (runNumber >= 104315 && runNumber <= 104321) {
484     return "4x4a*";
485   }
486   else if (runNumber >= 104792 && runNumber <= 104803) {
487     return "4x4b";
488   }
489   else if (runNumber >= 104824 && runNumber <= 104892) {
490     return "4x4c";
491   }
492   else if (runNumber == 105143 || runNumber == 105160) {
493     return "16x16a";
494   }
495   else if (runNumber >= 105256 && runNumber <= 105268) {
496     return "4x4c";
497   } else if (runNumber >= 114786 && runNumber <= 116684) {
498     return "Single_2b_1_1_1";
499   }
500   else {
501     AliError(Form("Unknown filling scheme (run %d)", runNumber));
502   }
503
504   return "Unknown";
505 }
506
507 Int_t AliPhysicsSelection::GetRatioBBBE(Int_t runNumber) {
508
509
510   if(fMC) return 1;
511
512   if (runNumber == 105143 || runNumber == 105160) {
513     return 8;
514   }else if (runNumber == 114786 || runNumber == 114798 ) {
515     return 1;
516   } else if (runNumber >= 114783 && runNumber <= 116684){
517     return 1;
518   }
519   else if (fComputeBG &&
520            !(runNumber >= 105256 && runNumber <= 105268) &&
521            !(runNumber >= 104065 && runNumber <= 104160) &&
522            !(runNumber >= 104315 && runNumber <= 104321) &&
523            !(runNumber >= 104792 && runNumber <= 104803) &&
524            !(runNumber >= 104824 && runNumber <= 104892)
525            ){     
526
527     AliError(Form("Unknown run %d, assuming ratio BE/EE = 2",runNumber));
528
529   }
530
531   return 2;
532 }
533
534
535 const char * AliPhysicsSelection::GetBXIDs(UInt_t runNumber, const char * trigger)  {
536
537   if (!fUseBXNumbers || fMC) return "";
538
539   if      (runNumber >= 104065 && runNumber <= 104160) {
540     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2128 #3019";
541     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #346 #3465";
542     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1234 #1680";
543     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
544     else AliError(Form("Unknown trigger: %s", trigger));
545   } 
546   else if (runNumber >= 104315 && runNumber <= 104321) {
547     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2000 #2891";
548     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #218 #3337";
549     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1106 #1552";
550     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
551     else AliError(Form("Unknown trigger: %s", trigger));
552   }
553   else if (runNumber >= 104792 && runNumber <= 104803) {
554     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2228 #3119";
555     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2554 #446";
556     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1334 #769";
557     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
558     else AliError(Form("Unknown trigger: %s", trigger));
559   }
560   else if (runNumber >= 104824 && runNumber <= 104892) {
561     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #3119 #769";
562     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2554 #446";
563     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1334 #2228";
564     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
565     else AliError(Form("Unknown trigger: %s", trigger));
566   }
567   else if (runNumber == 105143 || runNumber == 105160) {
568     fRatioBEEE = 8;
569     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #1337 #1418 #2228 #2309 #3119 #3200 #446 #527";
570     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #1580  #1742  #1904  #2066  #2630  #2792  #2954  #3362";
571     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "  #845  #1007  #1169   #1577 #3359 #3521 #119  #281 ";
572     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
573     else AliError(Form("Unknown trigger: %s", trigger));
574   }
575   else if (runNumber >= 105256 && runNumber <= 105268) {
576     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #3019 #669";
577     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2454 #346";
578     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1234 #2128";
579     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
580     else AliError(Form("Unknown trigger: %s", trigger));
581   } else if (runNumber >= 114786 && runNumber <= 116684) { // 7 TeV 2010, assume always the same filling scheme
582     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #346";
583     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2131";
584     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #3019";
585     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #1238";
586     else AliError(Form("Unknown trigger: %s", trigger));
587   }
588   else {
589     AliError(Form("Unknown run %d, using all BXs!",runNumber));
590   }
591
592   return "";
593 }
594     
595 Bool_t AliPhysicsSelection::Initialize(Int_t runNumber)
596 {
597   // initializes the object for the given run
598   // TODO having the run number here and parameters hardcoded is clearly temporary, a way needs to be found to have a CDB-like configuration also for analysis
599   
600   Bool_t oldStatus = TH1::AddDirectoryStatus();
601   TH1::AddDirectory(kFALSE);
602   
603   if(!fBin0CallBack) 
604     AliError("Bin0 Callback not set: will not fill the statistics for the bin 0");
605
606   if (fMC) {
607     // ovverride BX and bg options in case of MC
608     fComputeBG    = kFALSE;
609     fUseBXNumbers = kFALSE;
610   }
611
612   Int_t triggerScheme = GetTriggerScheme(runNumber);
613   if (!fUsingCustomClasses && fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
614     AliFatal("Processing several runs with different trigger schemes is not supported");
615   
616   if(fComputeBG && fCurrentRun != -1 && fCurrentRun != runNumber) 
617     AliFatal("Cannot process several runs because BG computation is requested");
618
619   if(fComputeBG && !fUseBXNumbers) 
620     AliFatal("Cannot compute BG id BX numbers are not used");
621   
622   if(fUseBXNumbers && fFillingScheme != "" && fFillingScheme != GetFillingScheme(runNumber))
623     AliFatal("Cannot process runs with different filling scheme if usage of BX numbers is requested");
624
625   fRatioBEEE          = GetRatioBBBE(runNumber);
626   fFillingScheme      = GetFillingScheme(runNumber);
627
628   if(fComputeBG) SetBIFactors(runNumber);
629
630   AliInfo(Form("Initializing for run %d", runNumber));
631   
632   // initialize first time?
633   if (fCurrentRun == -1)
634   {
635     if (fUsingCustomClasses) {
636       AliInfo("Using user-provided trigger classes");
637     } else {
638       switch (triggerScheme)
639       {
640       case 0:
641         fCollTrigClasses.Add(new TObjString(""));
642         break;
643         
644       case 1:
645         { // need a new scope to avoid cross-initialization errors
646
647           if (fUseMuonTriggers) {
648             // Muon trigger have the same BXIDs of the corresponding CINT triggers
649             fCollTrigClasses.Add(new TObjString(Form("%s%s ","+CMUS1B-ABCE-NOPF-MUON",  GetBXIDs(runNumber,"CINT1B-ABCE-NOPF-ALL"))));
650             fBGTrigClasses.Add  (new TObjString(Form("%s%s ","+CMUS1A-ABCE-NOPF-MUON",  GetBXIDs(runNumber,"CINT1A-ABCE-NOPF-ALL"))));
651             fBGTrigClasses.Add  (new TObjString(Form("%s%s ","+CMUS1C-ABCE-NOPF-MUON",  GetBXIDs(runNumber,"CINT1C-ABCE-NOPF-ALL"))));      
652             fBGTrigClasses.Add  (new TObjString(Form("%s%s ","+CMUS1-E-NOPF-MUON"    ,  GetBXIDs(runNumber,"CINT1-E-NOPF-ALL"))));
653           }
654           TObjString * cint1b = new TObjString(Form("%s%s","+CINT1B-ABCE-NOPF-ALL",  GetBXIDs(runNumber,"CINT1B-ABCE-NOPF-ALL")));
655           TObjString * cint1a = new TObjString(Form("%s%s","+CINT1A-ABCE-NOPF-ALL",  GetBXIDs(runNumber,"CINT1A-ABCE-NOPF-ALL")));
656           TObjString * cint1c = new TObjString(Form("%s%s","+CINT1C-ABCE-NOPF-ALL",  GetBXIDs(runNumber,"CINT1C-ABCE-NOPF-ALL")));
657           TObjString * cint1e = new TObjString(Form("%s%s","+CINT1-E-NOPF-ALL",      GetBXIDs(runNumber,"CINT1-E-NOPF-ALL"))    );
658           //
659           fCollTrigClasses.Add(cint1b);
660           fBGTrigClasses.Add(cint1a);
661           fBGTrigClasses.Add(cint1c);
662           fBGTrigClasses.Add(cint1e);
663
664         }
665         break;
666         
667       case 2:
668         fCollTrigClasses.Add(new TObjString("+CSMBB-ABCE-NOPF-ALL"));
669         fBGTrigClasses.Add(new TObjString("+CSMBA-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
670         fBGTrigClasses.Add(new TObjString("+CSMBC-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
671         break;
672         
673       default:
674         AliFatal(Form("Unsupported trigger scheme %d", triggerScheme));
675       }
676     }
677     
678     Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
679     
680     for (Int_t i=0; i<count; i++)
681     {
682       AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
683       triggerAnalysis->SetAnalyzeMC(fMC);
684       triggerAnalysis->EnableHistograms();
685       triggerAnalysis->SetSPDGFOThreshhold(1);
686       fTriggerAnalysis.Add(triggerAnalysis);
687     }
688       
689     // TODO: shall I really delete this?
690     if (fHistStatistics[0])
691       delete fHistStatistics[0];
692     if (fHistStatistics[1])
693       delete fHistStatistics[1];
694   
695     fHistStatistics[kStatIdxBin0] = BookHistStatistics("_Bin0");
696     fHistStatistics[kStatIdxAll]  = BookHistStatistics("");
697     
698     if (fHistBunchCrossing)
699       delete fHistBunchCrossing;
700   
701     fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5,  count, -0.5, -0.5 + count);
702
703     if (fHistTriggerPattern)
704       delete fHistTriggerPattern;
705     
706     const int ntrig=3;
707     Int_t n = 1;
708     const Int_t nbinTrig = TMath::Nint(TMath::Power(2,ntrig));
709
710     fHistTriggerPattern = new TH1F("fHistTriggerPattern", "Trigger pattern: FO + 2*v0A + 4*v0C", 
711                                    nbinTrig, -0.5, nbinTrig-0.5);    
712     fHistTriggerPattern->GetXaxis()->SetBinLabel(1,"NO TRIG");
713     fHistTriggerPattern->GetXaxis()->SetBinLabel(2,"FO");
714     fHistTriggerPattern->GetXaxis()->SetBinLabel(3,"v0A");
715     fHistTriggerPattern->GetXaxis()->SetBinLabel(4,"FO & v0A");
716     fHistTriggerPattern->GetXaxis()->SetBinLabel(5,"v0C");
717     fHistTriggerPattern->GetXaxis()->SetBinLabel(6,"FO & v0C");
718     fHistTriggerPattern->GetXaxis()->SetBinLabel(7,"v0A & v0C");
719     fHistTriggerPattern->GetXaxis()->SetBinLabel(8,"FO & v0A & v0C");
720
721   
722     n = 1;
723     for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
724     {
725       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
726       n++;
727     }
728     for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
729     {
730       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
731       n++;
732     }
733
734     
735
736   }
737     
738   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
739   for (Int_t i=0; i<count; i++)
740   {
741     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
742   
743     switch (runNumber)
744     {
745       case 104315:
746       case 104316:
747       case 104320:
748       case 104321:
749       case 104439:
750         triggerAnalysis->SetV0TimeOffset(7.5);
751         break;
752       default:
753         triggerAnalysis->SetV0TimeOffset(0);
754     }
755   }
756     
757   fCurrentRun = runNumber;
758   
759   TH1::AddDirectory(oldStatus);
760   
761   return kTRUE;
762 }
763
764 TH2F * AliPhysicsSelection::BookHistStatistics(const char * tag) {
765
766     // add 6 rows to count for the estimate of good, accidentals and
767     // BG and the ratio of BG and accidentals to total +ratio goot to
768     // first col + 2 for error on good.
769
770   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
771 #ifdef VERBOSE_STAT
772   Int_t extrarows = fComputeBG ? 8 : 0;
773 #else
774   Int_t extrarows = fComputeBG ? 3 : 0;
775 #endif
776   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);
777
778   h->GetXaxis()->SetBinLabel(kStatTriggerClass,  "Trigger class");
779   h->GetXaxis()->SetBinLabel(kStatHWTrig,        "Hardware trigger");
780   h->GetXaxis()->SetBinLabel(kStatFO1,           "FO >= 1");
781   h->GetXaxis()->SetBinLabel(kStatFO2,           "FO >= 2");
782   h->GetXaxis()->SetBinLabel(kStatV0A,           "V0A");
783   h->GetXaxis()->SetBinLabel(kStatV0C,           "V0C");
784   h->GetXaxis()->SetBinLabel(kStatFMD,           "FMD");
785   h->GetXaxis()->SetBinLabel(kStatSSD1,          "SSD >= 2");
786   h->GetXaxis()->SetBinLabel(kStatV0ABG,         "V0A BG");
787   h->GetXaxis()->SetBinLabel(kStatV0CBG,         "V0C BG");
788   h->GetXaxis()->SetBinLabel(kStatMB1,           "(FO >= 1 | V0A | V0C) & !V0 BG");
789   h->GetXaxis()->SetBinLabel(kStatMB1Prime,      "(FO >= 2 | (FO >= 1 & (V0A | V0C)) | (V0A &v0C) ) & !V0 BG");
790   h->GetXaxis()->SetBinLabel(kStatFO1AndV0,      "FO >= 1 & (V0A | V0C) & !V0 BG");
791   h->GetXaxis()->SetBinLabel(kStatV0,            "V0A & V0C & !V0 BG & !BG ID");
792   h->GetXaxis()->SetBinLabel(kStatOffline,       "Offline Trigger");
793   h->GetXaxis()->SetBinLabel(kStatAny2Hits,      "2 Hits & !V0 BG");
794   h->GetXaxis()->SetBinLabel(kStatBG,            "Background identification");
795   h->GetXaxis()->SetBinLabel(kStatAccepted,      "Accepted");
796
797   Int_t n = 1;
798   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
799     {
800       h->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
801       n++;
802     }
803   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
804     {
805       h->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
806       n++;
807     }
808
809   if(fComputeBG) {
810     h->GetYaxis()->SetBinLabel(n++, "BG (A+C)");
811     h->GetYaxis()->SetBinLabel(n++, "ACC");
812 #ifdef VERBOSE_STAT
813     h->GetYaxis()->SetBinLabel(n++, "BG (A+C) %  (rel. to CINT1B)");
814     h->GetYaxis()->SetBinLabel(n++, "ACC % (rel. to CINT1B)");
815     h->GetYaxis()->SetBinLabel(n++, "ERR GOOD %");
816     h->GetYaxis()->SetBinLabel(n++, "GOOD % (rel. to 1st col)");
817     h->GetYaxis()->SetBinLabel(n++, "ERR GOOD");
818 #endif
819     h->GetYaxis()->SetBinLabel(n++, "GOOD");
820   }
821
822   return h;
823 }
824
825 void AliPhysicsSelection::Print(Option_t* /* option */) const
826 {
827   // print the configuration
828   
829   Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
830   
831   Printf("Collision trigger classes:");
832   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
833     Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
834   
835   Printf("Background trigger classes:");
836   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
837     Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
838
839   AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
840   
841   if (triggerAnalysis)
842   {
843     if (triggerAnalysis->GetV0TimeOffset() > 0)
844       Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
845     
846     Printf("\nTotal available events:");
847     
848     triggerAnalysis->PrintTriggerClasses();
849   }
850   
851   if (fHistStatistics[kStatIdxAll] && fCollTrigClasses.GetEntries() > 0)
852   {
853     Printf("\nSelection statistics for first collision trigger (%s):", ((TObjString*) fCollTrigClasses.First())->String().Data());
854     
855     Printf("Total events with correct trigger class: %d", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(1, 1));
856     Printf("Selected collision candidates: %d", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(fHistStatistics[kStatIdxAll]->GetXaxis()->FindBin("Accepted"), 1));
857   }
858   
859   if (fHistBunchCrossing)
860   {
861     Printf("\nBunch crossing statistics:");
862     
863     for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
864     {
865       TString str;
866       str.Form("Trigger %s has accepted events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
867       
868       for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
869         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
870           str += Form("%d, ", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
871        
872       Printf("%s", str.Data());
873     }
874     
875     for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
876     {
877       Int_t count = 0;
878       for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
879       {
880         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
881           count++;
882       }
883       if (count > 1)
884         Printf("WARNING: Bunch crossing %d has more than one trigger class active. Check BPTX functioning for this run!", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
885     }
886   }
887
888   if (fUsingCustomClasses)        
889     Printf("WARNING: Using custom trigger classes!");
890   if (fSkipTriggerClassSelection) 
891     Printf("WARNING: Skipping trigger class selection!");
892   if (fSkipV0) 
893     Printf("WARNING: Ignoring V0 information in selection");
894   if(!fBin0CallBack) 
895     Printf("WARNING: Callback not set: will not fill the statistics for the bin 0");
896
897 }
898
899 Long64_t AliPhysicsSelection::Merge(TCollection* list)
900 {
901   // Merge a list of AliMultiplicityCorrection objects with this (needed for
902   // PROOF).
903   // Returns the number of merged objects (including this).
904
905   if (!list)
906     return 0;
907
908   if (list->IsEmpty())
909     return 1;
910
911   TIterator* iter = list->MakeIterator();
912   TObject* obj;
913   
914   // collections of all histograms
915   const Int_t nHists = 9;
916   TList collections[nHists];
917
918   Int_t count = 0;
919   while ((obj = iter->Next())) {
920
921     AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
922     if (entry == 0) 
923       continue;
924       
925     collections[0].Add(&(entry->fTriggerAnalysis));
926     if (entry->fHistStatistics[0])
927       collections[1].Add(entry->fHistStatistics[0]);
928     if (entry->fHistStatistics[1])
929       collections[2].Add(entry->fHistStatistics[1]);
930     if (entry->fHistBunchCrossing)
931       collections[3].Add(entry->fHistBunchCrossing);
932     if (entry->fHistTriggerPattern)
933       collections[4].Add(entry->fHistTriggerPattern);
934     if (entry->fBackgroundIdentification)
935       collections[5].Add(entry->fBackgroundIdentification);
936
937     count++;
938   }
939
940   fTriggerAnalysis.Merge(&collections[0]);
941   if (fHistStatistics[0])
942     fHistStatistics[0]->Merge(&collections[1]);
943   if (fHistStatistics[1])
944     fHistStatistics[1]->Merge(&collections[2]);
945   if (fHistBunchCrossing)
946     fHistBunchCrossing->Merge(&collections[3]);
947   if (fHistTriggerPattern)
948     fHistTriggerPattern->Merge(&collections[4]);
949   if (fBackgroundIdentification)
950     fBackgroundIdentification->Merge(&collections[5]);
951   
952   delete iter;
953
954   return count+1;
955 }
956
957 void AliPhysicsSelection::SaveHistograms(const char* folder) const
958 {
959   // write histograms to current directory
960   
961   if (!fHistStatistics[0] || !fHistStatistics[1])
962     return;
963     
964   if (folder)
965   {
966     gDirectory->mkdir(folder);
967     gDirectory->cd(folder);
968   }
969   
970
971   // Fill the last rows of fHistStatistics before saving
972   if (fComputeBG) {
973     Int_t triggerScheme = GetTriggerScheme(UInt_t(fCurrentRun));
974     if(triggerScheme != 1){
975       AliWarning("BG estimate only supported for trigger scheme \"1\" (CINT1 suite)");
976     } else if (fUseMuonTriggers) {
977       AliWarning("BG estimate with muon triggers to be implemented");
978     } else {
979       Int_t nHistStat = 2;
980       // TODO: get number of rows in a more flexible way
981       // 1. loop over all cols
982
983       for(Int_t iHistStat = 0; iHistStat < nHistStat; iHistStat++){
984     
985         Int_t ncol = fHistStatistics[iHistStat]->GetNbinsX();
986         Float_t good1 = 0;
987         for(Int_t icol = 1; icol <= ncol; icol++) {
988           Int_t cint1B = (Int_t) fHistStatistics[iHistStat]->GetBinContent(icol,1);     
989           Int_t cint1A = (Int_t) fHistStatistics[iHistStat]->GetBinContent(icol,2);     
990           Int_t cint1C = (Int_t) fHistStatistics[iHistStat]->GetBinContent(icol,3);     
991           Int_t cint1E = (Int_t) fHistStatistics[iHistStat]->GetBinContent(icol,4);      
992       
993           if (cint1B>0) {
994             Int_t acc  = fRatioBEEE*cint1E; 
995             Double_t acc_err = TMath::Sqrt(fRatioBEEE*fRatioBEEE*cint1E);
996             //      Int_t bg   = cint1A + cint1C - 2*acc;
997             Float_t bg   = fBIFactorA*(cint1A-acc) + fBIFactorC*(cint1C-acc) ;
998             Float_t good = Float_t(cint1B) - bg - acc;
999             if (icol ==1) good1 = good;
1000             //      Float_t errGood     = TMath::Sqrt(2*(cint1A+cint1C+cint1E));// Error on the number of goods assuming only bg fluctuates
1001             //      DeltaG^2 = B + FA^2 A + FC^2 C + Ratio^2 (FA+FC-1)^2 E.
1002             Float_t errGood     = TMath::Sqrt( cint1B + 
1003                                                fBIFactorA*fBIFactorA*cint1A +
1004                                                fBIFactorC*fBIFactorC*cint1C +
1005                                                fRatioBEEE * fRatioBEEE * 
1006                                                (fBIFactorA + fBIFactorC - 1)*(fBIFactorA + fBIFactorC - 1)*cint1E);
1007             Float_t errBG = TMath::Sqrt(fBIFactorA*fBIFactorA*cint1A+
1008                                         fBIFactorC*fBIFactorC*cint1C+
1009                                         fRatioBEEE*fRatioBEEE*(fBIFactorA+fBIFactorC)*(fBIFactorA+fBIFactorC)*cint1E);
1010         
1011         
1012             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowBG,bg);      
1013             fHistStatistics[iHistStat]->SetBinError  (icol,kStatRowBG,errBG);   
1014             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowAcc,acc);    
1015             fHistStatistics[iHistStat]->SetBinError  (icol,kStatRowAcc,acc_err);        
1016             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowGood,good);    
1017             fHistStatistics[iHistStat]->SetBinError  (icol,kStatRowGood,errGood);    
1018
1019 #ifdef VERBOSE_STAT
1020             //kStatRowBG=5,kStatRowAcc,kStatRowBGFrac,kStatRowAccFrac,kStatRowErrGoodFrac,kStatRowGoodFrac,kStatRowGood,kStatRowErrGood
1021             Float_t accFrac   = Float_t(acc) / cint1B  *100;
1022             Float_t errAccFrac= Float_t(acc_err) / cint1B  *100;
1023             Float_t bgFrac    = Float_t(bg)  / cint1B  *100;
1024             Float_t goodFrac  = Float_t(good)  / good1 *100;
1025             Float_t errGoodFrac = errGood/good1 * 100;
1026             Float_t errFracBG = bg > 0 ? TMath::Sqrt((errBG/bg)*(errBG/bg) + 1/cint1B)*bgFrac : 0;
1027             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowBGFrac,bgFrac);      
1028             fHistStatistics[iHistStat]->SetBinError  (icol,kStatRowBGFrac,errFracBG);   
1029             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowAccFrac,accFrac);    
1030             fHistStatistics[iHistStat]->SetBinError  (icol,kStatRowAccFrac,errAccFrac);    
1031             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowGoodFrac,goodFrac);    
1032             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowErrGoodFrac,errGoodFrac);    
1033             fHistStatistics[iHistStat]->SetBinContent(icol,kStatRowErrGood,errGood);    
1034 #endif
1035           }
1036         }
1037       }
1038     }  
1039   }
1040
1041   fHistStatistics[0]->Write();
1042   fHistStatistics[1]->Write();
1043   if(fHistBunchCrossing ) fHistBunchCrossing ->Write();
1044   if(fHistTriggerPattern) fHistTriggerPattern->Write();
1045   
1046   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
1047   for (Int_t i=0; i < count; i++)
1048   {
1049     TString triggerClass = "trigger_histograms_";
1050     if (i < fCollTrigClasses.GetEntries())
1051       triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
1052     else
1053       triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
1054   
1055     gDirectory->mkdir(triggerClass);
1056     gDirectory->cd(triggerClass);
1057   
1058     static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
1059     
1060     gDirectory->cd("..");
1061   }
1062  
1063   if (fBackgroundIdentification)
1064   {
1065     gDirectory->mkdir("background_identification");
1066     gDirectory->cd("background_identification");
1067       
1068     fBackgroundIdentification->GetOutput()->Write();
1069       
1070     gDirectory->cd("..");
1071   }
1072   
1073   if (folder)
1074     gDirectory->cd("..");
1075 }
1076
1077 void AliPhysicsSelection::SetBIFactors(Int_t run) {
1078
1079   switch(run) {
1080   case 104155:
1081     fBIFactorA = 0.961912722908;
1082     fBIFactorC = 1.04992336081;
1083     break;
1084   case 104157:
1085     fBIFactorA = 0.947312854998;
1086     fBIFactorC = 1.01599706417;
1087     break;
1088   case 104159:
1089     fBIFactorA = 0.93659320151;
1090     fBIFactorC = 0.98580804207;
1091     break;
1092   case 104160:
1093     fBIFactorA = 0.929664189926;
1094     fBIFactorC = 0.963467679851;
1095     break;
1096   case 104315:
1097     fBIFactorA = 1.08939104979;
1098     fBIFactorC = 0.931113921925;
1099     break;
1100   case 104316:
1101     fBIFactorA = 1.08351880974;
1102     fBIFactorC = 0.916068345845;
1103     break;
1104   case 104320:
1105     fBIFactorA = 1.07669281245;
1106     fBIFactorC = 0.876818744763;
1107     break;
1108   case 104321:
1109     fBIFactorA = 1.00971079602;
1110     fBIFactorC = 0.773781299076;
1111     break;
1112   case 104792:
1113     fBIFactorA = 0.787215863962;
1114     fBIFactorC = 0.778253173071;
1115     break;
1116   case 104793:
1117     fBIFactorA = 0.692211363661;
1118     fBIFactorC = 0.733152456667;
1119     break;
1120   case 104799:
1121     fBIFactorA = 1.04027825161;
1122     fBIFactorC = 1.00530825942;
1123     break;
1124   case 104800:
1125     fBIFactorA = 1.05309910671;
1126     fBIFactorC = 1.00376801855;
1127     break;
1128   case 104801:
1129     fBIFactorA = 1.0531231922;
1130     fBIFactorC = 0.992439666758;
1131     break;
1132   case 104802:
1133     fBIFactorA = 1.04191478134;
1134     fBIFactorC = 0.979368585208;
1135     break;
1136   case 104803:
1137     fBIFactorA = 1.03121314094;
1138     fBIFactorC = 0.973379962609;
1139     break;
1140   case 104824:
1141     fBIFactorA = 0.969945926722;
1142     fBIFactorC = 0.39549745806;
1143     break;
1144   case 104825:
1145     fBIFactorA = 0.968627213937;
1146     fBIFactorC = 0.310100412205;
1147     break;
1148   case 104841:
1149     fBIFactorA = 0.991601393212;
1150     fBIFactorC = 0.83762204722;
1151     break;
1152   case 104845:
1153     fBIFactorA = 0.98040863886;
1154     fBIFactorC = 0.694824205793;
1155     break;
1156   case 104867:
1157     fBIFactorA = 1.10646173412;
1158     fBIFactorC = 0.841407246916;
1159     break;
1160   case 104876:
1161     fBIFactorA = 1.12063452421;
1162     fBIFactorC = 0.78726542895;
1163     break;
1164   case 104890:
1165     fBIFactorA = 1.02346137453;
1166     fBIFactorC = 1.03355663595;
1167     break;
1168   case 104892:
1169     fBIFactorA = 1.05406025913;
1170     fBIFactorC = 1.00029166135;
1171     break;
1172   case 105143:
1173     fBIFactorA = 0.947343384349;
1174     fBIFactorC = 0.972637444408;
1175     break;
1176   case 105160:
1177     fBIFactorA = 0.908854622177;
1178     fBIFactorC = 0.958851103977;
1179     break; 
1180   case 105256:
1181     fBIFactorA = 0.810076150206;
1182     fBIFactorC = 0.884663561883;
1183     break;
1184   case 105257:
1185     fBIFactorA = 0.80974912303;
1186     fBIFactorC = 0.878859123479;
1187     break;
1188   case 105268:
1189     fBIFactorA = 0.809052110679;
1190     fBIFactorC = 0.87233890989;
1191     break;
1192   default:
1193     fBIFactorA = 1;
1194     fBIFactorC = 1;
1195   }
1196
1197
1198 }