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