]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/RESONANCES/macros/lego_train/old/AddRsnPairsKStar.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / lego_train / old / AddRsnPairsKStar.C
CommitLineData
900e06e4 1void AddRsnPairsKStar(AliAnalysisTaskSE *task,
2 Bool_t isMC,
3 Bool_t isMixing,
4 AliPID::EParticleType pType1,
5 Int_t listID1,
6 AliPID::EParticleType pType2,
7 Int_t listID2,
8 AliRsnCutSet *cutsEvent=0,
9 AliRsnCutSet *cutsPair=0,
10 TString suffix = "") {
11
12 Printf("id1=%d id2=%d",listID1,listID2);
13
14 // retrieve mass from PDG database
15 Int_t pdg = 313;
16 TDatabasePDG *db = TDatabasePDG::Instance();
17 TParticlePDG *part = db->GetParticle(pdg);
18 Double_t mass = part->Mass();
19
547e2d97 20 Bool_t valid;
21 Int_t isRsnMini = AliAnalysisManager::GetGlobalInt("rsnUseMiniPackage",valid);
22 if (isRsnMini) {
900e06e4 23 AddPairOutputMiniKStar(task,isMC,isMixing,pType1,listID1,pType2,listID2,pdg,mass,cutsPair,suffix);
24 } else {
25 // this function is common and it is located in RsnConfig.C
26 // as ouptup AddPairOutputKStar from this macro will be taken
27 AddPair(task,isMC,isMixing,pType1,listID1,pType2,listID2,pdg,mass,cutsEvent,cutsPair,suffix);
28 }
29}
30void AddPairOutputKStar(AliRsnLoopPair *pair)
31{
547e2d97 32 Bool_t valid;
33 Int_t isFullOutput = AliAnalysisManager::GetGlobalInt("rsnOutputFull",valid);
89ccdaab 34 Int_t isPP = AliAnalysisManager::GetGlobalInt("rsnIsPP",valid);
dc8dd82a 35
900e06e4 36 // axes
37 AliRsnValuePair *axisIM = new AliRsnValuePair("IM", AliRsnValuePair::kInvMass);
900e06e4 38 axisIM ->SetBins(900, 0.6, 1.5);
89ccdaab 39
40 AliRsnValuePair *axisPt = new AliRsnValuePair("PT", AliRsnValuePair::kPt);
900e06e4 41 axisPt ->SetBins(120, 0.0, 12.0);
42
89ccdaab 43 AliRsnValuePair *axisEta = new AliRsnValuePair("ETA", AliRsnValuePair::kEta);
44 axisEta ->SetBins(400, -0.5, 0.5);
45
46 AliRsnValueEvent *axisCentrality = 0;
47 if (!isPP) axisCentrality = new AliRsnValueEvent("MULTI",AliRsnValueEvent::kCentralityV0);
48
900e06e4 49 // output: 2D histogram of inv. mass vs. pt
50 AliRsnListOutput *outPair = 0;
547e2d97 51 if (!isFullOutput) {
900e06e4 52 outPair = new AliRsnListOutput("pair", AliRsnListOutput::kHistoDefault);
53 outPair->AddValue(axisIM);
54 } else {
55 outPair = new AliRsnListOutput("pair", AliRsnListOutput::kHistoSparse);
56 outPair->AddValue(axisIM);
57 outPair->AddValue(axisPt);
89ccdaab 58 outPair->AddValue(axisEta);
59 if (axisCentrality) axisCentrality->SetBins(20,0,100);
900e06e4 60 }
61 // add outputs to loop
62 pair->AddOutput(outPair);
63}
64
65void AddPairOutputMiniKStar(AliAnalysisTaskSE *task,Bool_t isMC,Bool_t isMixing, AliPID::EParticleType pType1,Int_t listID1, AliPID::EParticleType pType2,Int_t listID2, Int_t pdgMother,Double_t massMother, AliRsnCutSet *cutsPair=0,TString suffix = "") {
66
547e2d97 67 Bool_t valid;
68 Int_t isFullOutput = AliAnalysisManager::GetGlobalInt("rsnOutputFull",valid);
69 Int_t useMixing = AliAnalysisManager::GetGlobalInt("rsnUseMixing",valid);
89ccdaab 70 Int_t isPP = AliAnalysisManager::GetGlobalInt("rsnIsPP",valid);
dc8dd82a 71
900e06e4 72 AliRsnMiniAnalysisTask *taskRsnMini = (AliRsnMiniAnalysisTask *)task;
dc8dd82a 73
74 if (isPP) taskRsnMini->UseMultiplicity("QUALITY");
75 else {
76 taskRsnMini->UseCentrality("V0M");
77 Int_t multID = taskRsnMini->CreateValue(AliRsnMiniValue::kMult, kFALSE);
78 AliRsnMiniOutput *outMult = taskRsnMini->CreateOutput("eventMult", "HIST", "EVENT");
79 outMult->AddAxis(multID, 100, 0.0, 100.0);
80 Int_t paID = taskRsnMini->CreateValue(AliRsnMiniValue::kPlaneAngle, kFALSE);
81 AliRsnMiniOutput *outPa = taskRsnMini->CreateOutput("planeAngle", "HIST", "EVENT");
82 outPa->AddAxis(paID, 100, 0, TMath::Pi());
83 }
84
900e06e4 85 /* invariant mass */ Int_t imID = taskRsnMini->CreateValue(AliRsnMiniValue::kInvMass, kFALSE);
86 /* IM resolution */ Int_t resID = taskRsnMini->CreateValue(AliRsnMiniValue::kInvMassRes, kTRUE);
87 /* transv. momentum */ Int_t ptID = taskRsnMini->CreateValue(AliRsnMiniValue::kPt, kFALSE);
88 /* centrality */ Int_t centID = taskRsnMini->CreateValue(AliRsnMiniValue::kMult, kFALSE);
89ccdaab 89 /* eta */ Int_t etaID = taskRsnMini->CreateValue(AliRsnMiniValue::kEta, kFALSE);
dc8dd82a 90 /* rapidity */ Int_t yID = taskRsnMini->CreateValue(AliRsnMiniValue::kY, kFALSE);
91
92 Bool_t useRapidity = kTRUE;
900e06e4 93
dc8dd82a 94 Int_t nIM = 90; Double_t minIM = 0.6, maxIM = 1.5;
95 Int_t nRes = 200; Double_t minRes = -0.02, maxRes = 0.02;
96 Int_t nEta = 400; Double_t minEta = -0.5, maxEta = 0.5;
97 Int_t nY = 16; Double_t minY = -0.8, maxY = 0.8;
98 Int_t nPt = 120; Double_t minPt = 0.0, maxPt = 12.0;
99 Int_t nCent = 100; Double_t minCent = 0.0, maxCent = 100.0;
900e06e4 100 //
101 // -- Create all needed outputs -----------------------------------------------------------------
102 //
103
104 Int_t iCutK = listID1;
105 Int_t iCutPi = listID2;
106
547e2d97 107 // common definitions
108 TString outputType = "HIST";
109 if (isFullOutput) outputType = "SPARSE";
110
900e06e4 111 // use an array for more compact writing, which are different on mixing and charges
112 // [0] = unlike
113 // [1] = mixing
114 // [2] = like ++
115 // [3] = like --
dc8dd82a 116// const Int_t num = 12;
117 Bool_t use [12] = { 1 , 1 , useMixing, useMixing, 1 , 1 , isMC , isMC , isMC , isMC , isMC , isMC };
118 Bool_t useIM [12] = { 1 , 1 , 1 , 1 , 1 , 1 , 1 , 1 , 0 , 0 , 1 , 1 };
119 TString name [12] = {"Unlike1", "Unlike2", "Mixing1", "Mixing2", "LikePP", "LikeMM", "Trues1", "Trues2", "Res1" , "Res2" , "Mother1", "Mother2" };
120 TString comp [12] = {"PAIR" , "PAIR" , "MIX" , "MIX" , "PAIR" , "PAIR" , "TRUE" , "TRUE" , "TRUE" , "TRUE" , "MOTHER" , "MOTHER" };
121 Char_t charge1 [12] = {'+' , '-' , '+' , '-' , '+' , '-' , '+' , '-' , '+' , '-' , '+' , '-' };
122 Char_t charge2 [12] = {'-' , '+' , '-' , '+' , '+' , '-' , '-' , '+' , '-' , '+' , '-' , '+' };
123 Int_t cutID1 [12] = { iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK , iCutK };
124 Int_t cutID2 [12] = { iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi , iCutPi };
125
126 for (Int_t i = 0; i < 12; i++) {
900e06e4 127 if (!use[i]) continue;
128 // create output
547e2d97 129 AliRsnMiniOutput *out = taskRsnMini->CreateOutput(Form("%s_%s", suffix.Data(), name[i].Data()), outputType.Data(), comp[i].Data());
900e06e4 130 // selection settings
131 out->SetCutID(0, cutID1[i]);
132 out->SetCutID(1, cutID2[i]);
133 out->SetDaughter(0, AliRsnDaughter::kKaon);
134 out->SetDaughter(1, AliRsnDaughter::kPion);
135 out->SetCharge(0, charge1[i]);
136 out->SetCharge(1, charge2[i]);
137 out->SetMotherPDG(pdgMother);
138 out->SetMotherMass(massMother);
139 // pair cuts
140 if (cutsPair) out->SetPairCuts(cutsPair);
141 // axis X: invmass (or resolution)
142 if (useIM[i])
dc8dd82a 143 out->AddAxis(imID, nIM, minIM, maxIM);
900e06e4 144 else
dc8dd82a 145 out->AddAxis(resID, nRes, minRes, maxRes);
146
147 if (isFullOutput) {
148 // axis Y: transverse momentum
149 out->AddAxis(ptID, nPt, minPt, maxPt);
150 if (useRapidity) out->AddAxis(yID, nY, minY, maxY);
151 else out->AddAxis(etaID, nEta, minEta, maxEta);
152 // axis Z: centrality
153 if (!isPP) out->AddAxis(centID, nCent, minCent, maxCent);
154 }
155 }
156
157 // -- Create output for MC generated ------------------------------------------------------------
158 //
900e06e4 159
dc8dd82a 160 if (isMC) {
161 // create ouput
162 AliRsnMiniOutput *outMC = taskRsnMini->CreateOutput(Form("kstar_MCGen1%s", suffix.Data()), outputType.Data(), "MOTHER");
163 // selection settings
164 outMC->SetDaughter(0, AliRsnDaughter::kPion);
165 outMC->SetDaughter(1, AliRsnDaughter::kKaon);
166 outMC->SetMotherPDG(313);
167 outMC->SetMotherMass(massMother);
168 // pair cuts
169 if (cutsPair) outMC->SetPairCuts(cutsPair);
170 // axis X: invmass
171 outMC->AddAxis(imID, nIM, minIM, maxIM);
172 if (isFullOutput) {
173 // axis Y: transverse momentum
174 outMC->AddAxis(ptID, nPt, minPt, maxPt);
175 if (useRapidity) outMC->AddAxis(yID, nY, minY, maxY);
176 else outMC->AddAxis(etaID, nEta, minEta, maxEta);
177 // axis Z: centrality
178 if (!isPP) outMC->AddAxis(centID, nCent, minCent, maxCent);
179 }
180 }
181
182
183
184 if (isMC) {
185 // create ouput
186 AliRsnMiniOutput *outMC1 = taskRsnMini->CreateOutput(Form("phi_MCGen2%s", suffix.Data()), outputType.Data(), "MOTHER");
187 // selection settings
188 outMC1->SetDaughter(0, AliRsnDaughter::kKaon);
189 outMC1->SetDaughter(1, AliRsnDaughter::kPion);
190 outMC1->SetMotherPDG(-313);
191 outMC1->SetMotherMass(massMother);
192 // pair cuts
193 if (cutsPair) outMC1->SetPairCuts(cutsPair);
194 // axis X: invmass
195 outMC1->AddAxis(imID, nIM, minIM, maxIM);
547e2d97 196 if (isFullOutput) {
900e06e4 197 // axis Y: transverse momentum
dc8dd82a 198 outMC1->AddAxis(ptID, nPt, minPt, maxPt);
199 if (useRapidity) outMC1->AddAxis(yID, nY, minY, maxY);
200 else outMC1->AddAxis(etaID, nEta, minEta, maxEta);
900e06e4 201 // axis Z: centrality
dc8dd82a 202 if (!isPP) outMC1->AddAxis(centID, nCent, minCent, maxCent);
900e06e4 203 }
204 }
205
dc8dd82a 206
900e06e4 207}
dc8dd82a 208