]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/RESONANCES/macros/mini/ConfigTOFanalysisKStar.C
Added mother histo to Kstar config for TOF analysis + monitor histo binning changed
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / mini / ConfigTOFanalysisKStar.C
CommitLineData
72439a68 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 ConfigTOFanalysisKStar
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::kTOFpidKstarPbPb2010,
23 AliRsnCutSetDaughterParticle::ERsnDaughterCutSet cutKaCandidate = AliRsnCutSetDaughterParticle::kTOFpidKstarPbPb2010,
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
39 cutSetQ = new AliRsnCutSetDaughterParticle("cutQuality", AliRsnCutSetDaughterParticle::kQualityStd2010, AliPID::kPion, -1.0, aodFilterBit);
40 cutSetPi = new AliRsnCutSetDaughterParticle(Form("cutPionTOFPbPb2010_%2.1fsigma",nsigmaPi), cutPiCandidate, AliPID::kPion, nsigmaPi, aodFilterBit);
41 cutSetK = new AliRsnCutSetDaughterParticle(Form("cutKaonTOFPbPb2010_%2.1f2sigma",nsigmaKa), cutKaCandidate, AliPID::kKaon, nsigmaKa, aodFilterBit);
42
43 Int_t iCutQ = task->AddTrackCuts(cutSetQ);
44 Int_t iCutPi = task->AddTrackCuts(cutSetPi);
45 Int_t iCutK = task->AddTrackCuts(cutSetK);
46
47 if(enableMonitor){
48 Printf("======== Monitoring cut AliRsnCutSetDaughterParticle enabled");
da3e2cdb 49 gROOT->LoadMacro("$ALICE_ROOT/PWGLF/RESONANCES/macros/mini/AddMonitorOutput.C");
72439a68 50 AddMonitorOutput(isMC, cutSetQ->GetMonitorOutput());
51 AddMonitorOutput(isMC, cutSetPi->GetMonitorOutput());
52 AddMonitorOutput(isMC, cutSetK->GetMonitorOutput());
53 }
54
55 // -- Values ------------------------------------------------------------------------------------
56 /* invariant mass */ Int_t imID = task->CreateValue(AliRsnMiniValue::kInvMass, kFALSE);
57 /* IM resolution */ Int_t resID = task->CreateValue(AliRsnMiniValue::kInvMassRes, kTRUE);
58 /* transv. momentum */ Int_t ptID = task->CreateValue(AliRsnMiniValue::kPt, kFALSE);
59 /* centrality */ Int_t centID = task->CreateValue(AliRsnMiniValue::kMult, kFALSE);
60 /* pseudorapidity */ Int_t etaID = task->CreateValue(AliRsnMiniValue::kEta, kFALSE);
61 /* rapidity */ Int_t yID = task->CreateValue(AliRsnMiniValue::kY, kFALSE);
62
63 // -- Create all needed outputs -----------------------------------------------------------------
64 // use an array for more compact writing, which are different on mixing and charges
65 // [0] = unlike
66 // [1] = mixing
67 // [2] = like ++
68 // [3] = like --
69 Bool_t use [10] = { !IsMcTrueOnly, !IsMcTrueOnly, !IsMcTrueOnly, !IsMcTrueOnly , !IsMcTrueOnly, !IsMcTrueOnly, isMC , isMC , isMC , isMC };
70 Bool_t useIM [10] = { 1 , 1 , 1 , 1 , 1 , 1 , 1 , 1 , 0 , 0 };
71 TString name [10] = {"UnlikePM", "UnlikeMP", "MixingPM", "MixingMP", "LikePP", "LikeMM", "TruesPM", "TruesMP", "ResPM" , "ResMP" };
72 TString comp [10] = {"PAIR" , "PAIR" , "MIX" , "MIX" , "PAIR" , "PAIR" , "TRUE" , "TRUE" , "TRUE" , "TRUE" };
73 //TString output [10] = {"HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" , "HIST" };
74 TString output [10] = {"SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" , "SPARSE" };
75 Char_t charge1 [10] = {'+' , '-' , '+' , '-' , '+' , '-' , '+' , '-' , '+' , '-' };
76 Char_t charge2 [10] = {'-' , '+' , '-' , '+' , '+' , '-' , '-' , '+' , '-' , '+' };
77 Int_t cutID1 [10] = { iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK };
78 Int_t cutID2 [10] = { iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi };
79
80 for (Int_t i = 0; i < 10; i++) {
81 if (!use[i]) continue;
82 AliRsnMiniOutput *out = task->CreateOutput(Form("kstar_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
83 out->SetCutID(0, cutID1[i]);
84 out->SetCutID(1, cutID2[i]);
85 out->SetDaughter(0, AliRsnDaughter::kKaon);
86 out->SetDaughter(1, AliRsnDaughter::kPion);
87 out->SetCharge(0, charge1[i]);
88 out->SetCharge(1, charge2[i]);
89 out->SetMotherPDG(313);
90 out->SetMotherMass(0.89594);
91 out->SetPairCuts(cutsPair);
92
93 // axis X: invmass (or resolution)
94 if (useIM[i])
95 out->AddAxis(imID, 90, 0.6, 1.5);
96 else
97 out->AddAxis(resID, 200, -0.02, 0.02);
98
99 // axis Y: transverse momentum
100 out->AddAxis(ptID, 100, 0.0, 10.0);
101
102 // axis Z: centrality-multiplicity
103 if (!isPP)
104 out->AddAxis(centID, 100, 0.0, 100.0);
105 else
106 out->AddAxis(centID, 400, 0.0, 400.0);
107
108 // axis W: pseudorapidity
109 // out->AddAxis(etaID, 20, -1.0, 1.0);
110 // axis J: rapidity
111 // out->AddAxis(yID, 10, -0.5, 0.5);
112
113 }
58a1c9a4 114
115 if (isMC){
116 // create output
117 AliRsnMiniOutput *outm = task->CreateOutput(Form("kstar_Mother%s", suffix), "SPARSE", "MOTHER");
118 // selection settings
119 outm->SetDaughter(0, AliRsnDaughter::kKaon);
120 outm->SetDaughter(1, AliRsnDaughter::kPion);
121 outm->SetMotherPDG(313);
122 outm->SetMotherMass(0.89594);
123 // pair cuts
124 outm->SetPairCuts(cutsPair);
125 // binnings
126 outm->AddAxis(imID, 90, 0.6, 1.5);
127 outm->AddAxis(ptID, 100, 0.0, 10.0);
128
129 if (!isPP)
130 outm->AddAxis(centID, 100, 0.0, 100.0);
131 else
132 outm->AddAxis(centID, 400, 0.0, 400.0);
133 }
72439a68 134 return kTRUE;
135}