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