]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/AliD0toKpiPlots.C
Added missing getters (Andrea)
[u/mrichter/AliRoot.git] / PWG3 / AliD0toKpiPlots.C
1 void AliD0toKpiPlots(const Char_t *inName="AliD0toKpi.root",
2                      const Char_t *outName="D0histograms.root") {
3   //--------------------------------------------------------------------------
4   // This macro histograms many variables of D0->Kpi candidates
5   //
6   //     Andrea Dainese, andrea.dainese@lnl.infn.it
7   //--------------------------------------------------------------------------
8
9   gSystem->Load("libANALYSIS");
10   gSystem->Load("libAOD.so");
11   gSystem->Load("libPWG3base.so");
12
13   // set of cuts
14   Double_t D0Cuts[9] = {0.1,         // mass [GeV]
15                         1000000.,          // dca [micron]
16                         1.1,           // cosThetaStar
17                         0.,           // pT K [GeV/c]
18                         0.,           // pT Pi [GeV/c]    
19                         100000.,       // d0K upper [micron]
20                         100000.,       // d0Pi upper [micron]
21                         10000000000.,            // d0d0 [micron^2]
22                         -1.1};          // cosThetaPointing
23
24   // number of events (for normalization)
25   Bool_t normalize = kFALSE;
26   Double_t  events = 1.;
27
28
29
30   // define histograms
31   TH1F *hptK = new TH1F("hptK","\"K\" p_{t} distribution",50,0,10);
32   hptK->SetXTitle("p_{t} [GeV]");
33
34   TH1F *hptPi = new TH1F("hptPi","\"#pi\" p_{t} distribution",50,0,10);
35   hptPi->SetXTitle("p_{t} [GeV]");
36
37   TH1F *hDCA = new TH1F("hDCA","DCA",50,0,1000);
38   hDCA->SetXTitle("dca [#mu m]");
39
40   TH1F *hptD0 = new TH1F("hptD0","D^{0} p_{t} distribution",40,0,40);
41   hptD0->SetXTitle("p_{t} [GeV]");
42
43   TH1F *hyD0 = new TH1F("hyD0","D^{0} rapidity distribution",50,-2,2);
44   hyD0->SetXTitle("y");
45
46   TH1F *hCPtaD0 = new TH1F("hCPtaD0","cosine of pointing angle distribution",100,-1,1);
47   hCPtaD0->SetXTitle("cos #theta_{point}");
48
49   TH1F *hCPtaXY = new TH1F("hCPtaXY","cosine of pointing angle in (x,y) plane",100,-1,1);
50   hCPtaXY->SetXTitle("cos #theta_{point}");
51
52   TH1F *hCts = new TH1F("hCts","cosine of decay angle",50,-1.2,1.2);
53   hCts->SetXTitle("cos #theta^{*}");
54
55   TH2F *hCtsVsPtK = new TH2F("hCtsVsPtK","cosine of decay angle VS \"K\" p_{t}",50,0,5,50,-1,1);
56   hCtsVsPtK->SetYTitle("cos #theta^{*}");
57   hCtsVsPtK->SetXTitle("p_{t} [GeV]");
58
59   TH1F *hd0d0 = new TH1F("hd0d0","Product of the impact parameters",100,-100000,100000);
60   hd0d0->SetXTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
61
62   TH1F *hd0K = new TH1F("hd0K","Impact parameter of \"K\"",100,-5000,5000);
63   hd0K->SetXTitle("d_{0}^{K} [#mu m]");
64
65   TH1F *hd0Pi = new TH1F("hd0Pi","Impact parameter of \"#pi\"",100,-5000,5000);
66   hd0Pi->SetXTitle("d_{0}^{#pi} [#mu m]");
67
68   TH2F *hCPtaVsd0d0 = new TH2F("hCPtaVsd0d0","cos #theta_{point} vs d_{0}^{K} #times d_{0}^{#pi}",100,-100000,100000,100,-1,1);
69   hCPtaVsd0d0->SetXTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
70   hCPtaVsd0d0->SetYTitle("cos #theta_{point}");
71
72   TH2F *hCPtaVsd0d0zoom = new TH2F("hCPtaVsd0d0zoom","cos #theta_{point} vs d_{0}^{K} #times d_{0}^{#pi}",100,-100000,0,100,.9,1);
73   hCPtaVsd0d0zoom->SetXTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
74   hCPtaVsd0d0zoom->SetYTitle("cos #theta_{point}");
75
76   TH2F *hd0d0VSptD0 = new TH2F("hd0d0VSptD0","d_{0}^{K} #times d_{0}^{#pi} VS D^{0} p_{t}",50,0,25,100,-120000,120000);
77   hd0d0VSptD0->SetYTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
78   hd0d0VSptD0->SetXTitle("D^{0} p_{t} [GeV]");
79
80   TH1F *hMass = new TH1F("hMass","Invariant mass distribution",50,1.765,1.965);
81   hMass->SetXTitle("M[K,#pi] [GeV]");
82
83   TH2F *hArm = new TH2F("hArm","Armenteros plot",50,-2,2,50,0,1);
84   hArm->SetXTitle("#alpha");
85   hArm->SetYTitle("q_{t}");
86
87   // open input file and get tree
88   TFile *inFile = TFile::Open(inName);
89
90   TTree *treeD0 = (TTree*)inFile->Get("TreeD0");
91   AliD0toKpi *D = 0; 
92   treeD0->SetBranchAddress("D0toKpi",&D);
93   Int_t entries = (Int_t)treeD0->GetEntries();
94
95   printf("+++\n+++ Number of D0 in tree:  %d\n+++\n",entries);
96
97   Double_t MD0,MD0bar,ctsD0,ctsD0bar,ctsPiD0,ctsPiD0bar;
98   Double_t WgtD0,WgtD0bar;
99   Double_t sampleABC=0.;
100   Int_t okD0=0,okD0bar=0;
101   Int_t nSel = 0;
102   Int_t ptbin;
103
104   // loop on D0
105   for(Int_t i=0; i<entries; i++) {
106     if(i%10000==0) printf(" candidate %d of %d\n",i,entries);
107
108     // get event from tree
109     treeD0->GetEvent(i);
110     //--- select the PID strategy & compute weights
111     //    D->ApplyPID("TOFparam_PbPb");
112     //    D->ComputeWgts();
113     // get weights for the three samples A+B+C 
114     //    D->GetWgts(WgtD0,WgtD0bar,"ABC");
115     WgtD0 = 1.; WgtD0bar = 1.;
116
117     // normalize to 1 event
118     if(normalize) { WgtD0 /= events; WgtD0bar /= events; }
119
120     // check if candidate passes selection (as D0 or D0bar)
121     D->Select(D0Cuts,okD0,okD0bar);
122
123     // set weights to 0 if the candidate doesn't pass selection
124     if(!okD0)    WgtD0=0.; 
125     if(!okD0bar) WgtD0bar=0.;
126     if(okD0 || okD0bar) nSel++;
127
128     // count selected candidates
129     sampleABC += WgtD0 + WgtD0bar;
130
131     // inv mass and cosThetaStar
132     D->InvMass(MD0,MD0bar);
133     D->CosThetaStar(ctsD0,ctsD0bar);
134     
135     // fill histograms
136     hptK->Fill(D->PtChild(1),WgtD0);
137     hptK->Fill(D->PtChild(0),WgtD0bar);
138     hptPi->Fill(D->PtChild(0),WgtD0);
139     hptPi->Fill(D->PtChild(1),WgtD0bar);
140     hd0K->Fill(D->Getd0Child(1),WgtD0);
141     hd0K->Fill(D->Getd0Child(0),WgtD0bar);
142     hd0Pi->Fill(D->Getd0Child(0),WgtD0);
143     hd0Pi->Fill(D->Getd0Child(1),WgtD0bar);   
144     hMass->Fill(MD0,WgtD0);
145     hMass->Fill(MD0bar,WgtD0bar);
146     hCts->Fill(ctsD0,WgtD0);
147     hCts->Fill(ctsD0bar,WgtD0bar);
148     hCtsVsPtK->Fill(D->PtChild(1),ctsD0,WgtD0);
149     hCtsVsPtK->Fill(D->PtChild(0),ctsD0bar,WgtD0bar);
150     hDCA->Fill(D->GetDCA(),WgtD0+WgtD0bar);
151     hptD0->Fill(D->Pt(),WgtD0+WgtD0bar);
152     hyD0->Fill(D->Rapidity(),WgtD0+WgtD0bar);
153     hd0d0->Fill(D->ProdImpParams(),WgtD0+WgtD0bar);
154     hCPtaD0->Fill(D->CosPointing(),WgtD0+WgtD0bar);
155     hCPtaXY->Fill(D->CosPointingXY(),WgtD0+WgtD0bar);
156     hCPtaVsd0d0->Fill(D->ProdImpParams(),D->CosPointing(),WgtD0+WgtD0bar);
157     hd0d0VSptD0->Fill(D->Pt(),D->ProdImpParams(),WgtD0+WgtD0bar);
158     hCPtaVsd0d0zoom->Fill(D->ProdImpParams(),D->CosPointing(),WgtD0+WgtD0bar);
159     hArm->Fill(D->Alpha(),D->Qt(),WgtD0+WgtD0bar);
160     
161    
162   } // end loop on D0 candidates
163
164   inFile->Close();
165
166   printf("\n\n --- Total number of candidates passing selection: %d\n\n --- Sum of weights sample A+B+C: %f\n\n",nSel,sampleABC); 
167   
168   // draw histograms
169   TCanvas *c1 = new TCanvas("c1","pt K & pi",0,0,700,700);
170   c1->SetLogy(); 
171   hptK->Draw();
172   hptPi->Draw("same");
173
174   TCanvas *c2 = new TCanvas("c2","pt D0",0,0,700,700);
175   c2->SetLogy(); 
176   hptD0->Draw();
177
178   TCanvas *c3 = new TCanvas("c3","rapidity D0",0,0,700,700);
179   hyD0->Draw(); 
180
181   TCanvas *c4 = new TCanvas("c4","pointing angle",0,0,700,700);
182   hCPtaD0->Draw();
183
184   TCanvas *c5 = new TCanvas("c5","d0 x d0",0,0,700,700);
185   c5->SetLogy();
186   hd0d0->Draw();
187
188   TCanvas *c6 = new TCanvas("c6","pointing angle VS d0d0",0,0,700,700);
189   c6->SetLogz();
190   hCPtaVsd0d0->Draw("box");
191
192   TCanvas *c7 = new TCanvas("c7","mass",0,0,700,700);
193   hMass->Draw();
194
195   TCanvas *c8 = new TCanvas("c8","armenteros",0,0,700,700);
196   hArm->Draw("box");
197
198   TCanvas *c9 = new TCanvas("c9","decay angle",0,0,700,700); 
199   hCts->Draw();
200
201   TCanvas *c10 = new TCanvas("c10","dca",0,0,700,700);
202   c10->SetLogy();
203   hDCA->Draw();
204
205   TCanvas *c11 = new TCanvas("c11","d0 K & pi",0,0,700,700);
206   c11->SetLogy();
207   hd0K->Draw();
208   hd0Pi->Draw("same");
209
210   // write all histograms to file
211   TFile *outFile = new TFile(outName,"recreate");
212   hMass->Write();
213   hDCA->Write(); 
214   hCts->Write();
215   hCtsVsPtK->Write();
216   hArm->Write();
217   hCPtaVsd0d0->Write();
218   hd0d0VSptD0->Write();
219   hCPtaVsd0d0zoom->Write();
220   hd0d0->Write();
221   hCPtaD0->Write();
222   hCPtaXY->Write();
223   hptK->Write();
224   hptPi->Write();
225   hptD0->Write();
226   hyD0->Write();
227   hd0K->Write();
228   hd0Pi->Write();
229   outFile->Close();
230   
231   return;
232 }
233
234
235
236