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