]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/RESONANCES/macros/mini/ConfigTPCanalysisKStar.C
Fixed variable type for trigger mask + removed dependency from unused class
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / mini / ConfigTPCanalysisKStar.C
CommitLineData
6184d9a9 1/***************************************************************************
2 fbellini@cern.ch - last modified on 06/08/2012
3
4// *** Configuration script for K*, anti-K* analysis with 2010 PbPb runs ***
5//
6// A configuration script for RSN package needs to define the followings:
7//
8// (1) decay tree of each resonance to be studied, which is needed to select
9// true pairs and to assign the right mass to all candidate daughters
10// (2) cuts at all levels: single daughters, tracks, events
11// (3) output objects: histograms or trees
12****************************************************************************/
13
14Bool_t ConfigTPCanalysisKStar
15(
16 AliRsnMiniAnalysisTask *task,
17 Bool_t isMC,
18 Bool_t isPP,
19 const char *suffix,
20 AliRsnCutSet *cutsPair,
21 Int_t aodFilterBit = 5,
22 AliRsnCutSetDaughterParticle::ERsnDaughterCutSet cutPiCandidate = AliRsnCutSetDaughterParticle::kFastTPCpidNsigma,
23 AliRsnCutSetDaughterParticle::ERsnDaughterCutSet cutKaCandidate = AliRsnCutSetDaughterParticle::kFastTPCpidNsigma,
24 Float_t nsigmaPi = 2.0,
25 Float_t nsigmaKa = 2.0,
26 Bool_t enableMonitor = kTRUE,
27 Bool_t IsMcTrueOnly = kFALSE,
28 Int_t aodN = 0
29)
30{
31 // manage suffix
32 if (strlen(suffix) > 0) suffix = Form("_%s", suffix);
33
34 // set daughter cuts
35 AliRsnCutSetDaughterParticle * cutSetQ;
36 AliRsnCutSetDaughterParticle * cutSetPi;
37 AliRsnCutSetDaughterParticle * cutSetK;
38
7fc41557 39 //2010 cuts
40 //cutSetQ = new AliRsnCutSetDaughterParticle("cutQuality", AliRsnCutSetDaughterParticle::kQualityStd2010, AliPID::kPion, -1.0, aodFilterBit);
41 //cutSetPi = new AliRsnCutSetDaughterParticle(Form("cutPionTPCPbPb2010_%2.1fsigma",nsigmaPi), cutPiCandidate, AliPID::kPion, nsigmaPi, aodFilterBit);
42 //cutSetK = new AliRsnCutSetDaughterParticle(Form("cutKaonTPCPbPb2010_%2.1f2sigma",nsigmaKa), cutKaCandidate, AliPID::kKaon, nsigmaKa, aodFilterBit);
43
44 //2011 High-pT cuts
45 cutSetQ = new AliRsnCutSetDaughterParticle(Form("cutQ_bit%i",aodFilterBit), AliRsnCutSetDaughterParticle::kQualityStd2011, AliPID::kPion, -1.0, aodFilterBit, kTRUE);
46 cutSetQ->SetUse2011StdQualityCutsHighPt(kTRUE);
47
48 cutSetPi = new AliRsnCutSetDaughterParticle(Form("cutPi%i_%2.1fsigma",cutPiCandidate, nsigmaPi), cutPiCandidate, AliPID::kPion, nsigmaPi, aodFilterBit,kTRUE);
49 cutSetPi->SetUse2011StdQualityCutsHighPt(kTRUE);
50
51 cutSetK = new AliRsnCutSetDaughterParticle(Form("cutK%i_%2.1fsigma",cutPiCandidate, nsigmaKa), cutKaCandidate, AliPID::kKaon, nsigmaKa, aodFilterBit,kTRUE);
52 cutSetK->SetUse2011StdQualityCutsHighPt(kTRUE);
53
6184d9a9 54
55 Int_t iCutQ = task->AddTrackCuts(cutSetQ);
56 Int_t iCutPi = task->AddTrackCuts(cutSetPi);
57 Int_t iCutK = task->AddTrackCuts(cutSetK);
7fc41557 58
6184d9a9 59
60 if(enableMonitor){
61 Printf("======== Monitoring cut AliRsnCutSetDaughterParticle enabled");
62 gROOT->LoadMacro("$ALICE_ROOT/PWGLF/RESONANCES/macros/mini/AddMonitorOutput.C");
63 AddMonitorOutput(isMC, cutSetQ->GetMonitorOutput());
64 AddMonitorOutput(isMC, cutSetPi->GetMonitorOutput());
65 AddMonitorOutput(isMC, cutSetK->GetMonitorOutput());
66 }
67
68 // -- Values ------------------------------------------------------------------------------------
69 /* invariant mass */ Int_t imID = task->CreateValue(AliRsnMiniValue::kInvMass, kFALSE);
70 /* IM resolution */ Int_t resID = task->CreateValue(AliRsnMiniValue::kInvMassRes, kTRUE);
71 /* transv. momentum */ Int_t ptID = task->CreateValue(AliRsnMiniValue::kPt, kFALSE);
72 /* centrality */ Int_t centID = task->CreateValue(AliRsnMiniValue::kMult, kFALSE);
73 /* pseudorapidity */ Int_t etaID = task->CreateValue(AliRsnMiniValue::kEta, kFALSE);
74 /* rapidity */ Int_t yID = task->CreateValue(AliRsnMiniValue::kY, kFALSE);
75
76 // -- Create all needed outputs -----------------------------------------------------------------
77 // use an array for more compact writing, which are different on mixing and charges
78 // [0] = unlike
79 // [1] = mixing
80 // [2] = like ++
81 // [3] = like --
82 Bool_t use [10] = { !IsMcTrueOnly, !IsMcTrueOnly, !IsMcTrueOnly, !IsMcTrueOnly , !IsMcTrueOnly, !IsMcTrueOnly, isMC , isMC , isMC , isMC };
83 Bool_t useIM [10] = { 1 , 1 , 1 , 1 , 1 , 1 , 1 , 1 , 0 , 0 };
84 TString name [10] = {"UnlikePM", "UnlikeMP", "MixingPM", "MixingMP", "LikePP", "LikeMM", "TruesPM", "TruesMP", "ResPM" , "ResMP" };
85 TString comp [10] = {"PAIR" , "PAIR" , "MIX" , "MIX" , "PAIR" , "PAIR" , "TRUE" , "TRUE" , "TRUE" , "TRUE" };
86 //TString output [10] = {"HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" };
87 TString output [10] = {"SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" };
88 Char_t charge1 [10] = {'+' , '-' , '+' , '-' , '+' , '-' , '+' , '-' , '+' , '-' };
89 Char_t charge2 [10] = {'-' , '+' , '-' , '+' , '+' , '-' , '-' , '+' , '-' , '+' };
90 Int_t cutID1 [10] = { iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK };
91 Int_t cutID2 [10] = { iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi };
92
93 for (Int_t i = 0; i < 10; i++) {
94 if (!use[i]) continue;
95 AliRsnMiniOutput *out = task->CreateOutput(Form("kstar_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
96 out->SetCutID(0, cutID1[i]);
97 out->SetCutID(1, cutID2[i]);
98 out->SetDaughter(0, AliRsnDaughter::kKaon);
99 out->SetDaughter(1, AliRsnDaughter::kPion);
100 out->SetCharge(0, charge1[i]);
101 out->SetCharge(1, charge2[i]);
102 out->SetMotherPDG(313);
103 out->SetMotherMass(0.89594);
104 out->SetPairCuts(cutsPair);
105
106 // axis X: invmass (or resolution)
107 if (useIM[i])
108 out->AddAxis(imID, 90, 0.6, 1.5);
109 else
110 out->AddAxis(resID, 200, -0.02, 0.02);
111
112 // axis Y: transverse momentum
7fc41557 113 out->AddAxis(ptID, 300, 0.0, 30.0);
6184d9a9 114
115 // axis Z: centrality-multiplicity
116 if (!isPP)
117 out->AddAxis(centID, 100, 0.0, 100.0);
118 else
119 out->AddAxis(centID, 400, 0.0, 400.0);
120
121 // axis W: pseudorapidity
122 // out->AddAxis(etaID, 20, -1.0, 1.0);
123 // axis J: rapidity
124 // out->AddAxis(yID, 10, -0.5, 0.5);
125
126 }
127 return kTRUE;
128}