]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/KINK/AliResonanceKinkLikeSign.cxx
Adding missing include
[u/mrichter/AliRoot.git] / PWG2 / KINK / AliResonanceKinkLikeSign.cxx
CommitLineData
10eaad41 1/**************************************************************************
2 * Author: Paraskevi Ganoti, University of Athens (pganoti@phys.uoa.gr) *
3 * Contributors are mentioned in the code where appropriate. *
4 * *
5 * Permission to use, copy, modify and distribute this software and its *
6 * documentation strictly for non-commercial purposes is hereby granted *
7 * without fee, provided that the above copyright notice appears in all *
8 * copies and that both the copyright notice and this permission notice *
9 * appear in the supporting documentation. The authors make no claims *
10 * about the suitability of this software for any purpose. It is *
11 * provided "as is" without express or implied warranty. *
12 **************************************************************************/
13
14//----------------------------------------------------------------------------------------------------------------
15// class AliResonanceKinkLikeSign
16// Example of an analysis task for producing a like-sign background for resonances having at least one
17// kaon-kink in their decay products.
18// Background is calculated from a negative kaon kink and a negative track.
19//-----------------------------------------------------------------------------------------------------------------
20
21#include "TChain.h"
22#include "TTree.h"
23#include "TVector3.h"
24#include "TF1.h"
25#include "TH1D.h"
26#include "TLorentzVector.h"
27#include "AliAnalysisTaskSE.h"
28#include "AliAnalysisManager.h"
29
30#include "AliESDInputHandler.h"
31
32#include "AliResonanceKinkLikeSign.h"
33#include "AliESDkink.h"
34
35ClassImp(AliResonanceKinkLikeSign)
36
37//________________________________________________________________________
38AliResonanceKinkLikeSign::AliResonanceKinkLikeSign()
39 : AliAnalysisTaskSE(), fESD(0), fListOfHistos(0), f1(0), f2(0), fNegKaonLikeSign(0)
40{
41 // Constructor
42}
43
44//________________________________________________________________________
45AliResonanceKinkLikeSign::AliResonanceKinkLikeSign(const char *name)
46 : AliAnalysisTaskSE(name), fESD(0), fListOfHistos(0), f1(0), f2(0), fNegKaonLikeSign(0)
47
48{
49 // Constructor
50
51 // Define input and output slots here
52 // Input slot #0 works with a TChain
53 DefineInput(0, TChain::Class());
54 DefineOutput(1, TList::Class());
55}
56
57//________________________________________________________________________
58void AliResonanceKinkLikeSign::ConnectInputData(Option_t *)
59{
60 // Connect ESD or AOD here
61 // Called once
10eaad41 62
63 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
64 if (!tree) {
65 Printf("ERROR: Could not read chain from input slot 0");
66 } else {
67 // Disable all branches, we want to process only MC
68 tree->SetBranchStatus("*",kTRUE );
69 tree->SetBranchStatus("*Calo*", kFALSE);
70
71 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
72
73 if (!esdH) {
74 Printf("ERROR: Could not get ESDInputHandler");
75 } else
76 fESD = esdH->GetEvent();
77 }
78}
79
80//________________________________________________________________________
81void AliResonanceKinkLikeSign::CreateOutputObjects()
82{
83 // Create histograms
84 // Called once
85
86 f1=new TF1("f1","((atan([0]*[1]*(1.0/(sqrt((x^2)*(1.0-([1]^2))-([0]^2)*([1]^2))))))*180.)/[2]",1.1,10.0);
87 f1->SetParameter(0,0.493677);
88 f1->SetParameter(1,0.9127037);
89 f1->SetParameter(2,TMath::Pi());
90
91 f2=new TF1("f2","((atan([0]*[1]*(1.0/(sqrt((x^2)*(1.0-([1]^2))-([0]^2)*([1]^2))))))*180.)/[2]",0.1,10.0);
92 f2->SetParameter(0,0.13957018);
93 f2->SetParameter(1,0.2731374);
94 f2->SetParameter(2,TMath::Pi());
95
96 //OpenFile(1); //uncomment for proof analysis
97
98 // for K*0(892)
99 fNegKaonLikeSign=new TH1D("fNegKaonLikeSign"," ", 60,0.6,1.2);
100
101 // for phi(1020)
102 // fNegKaonLikeSign=new TH1D("fNegKaonLikeSign"," ", 70,0.99,1.088);
103
104
105 fListOfHistos=new TList();
106 fListOfHistos->Add(fNegKaonLikeSign);
107
108}
109
110//________________________________________________________________________
111void AliResonanceKinkLikeSign::Exec(Option_t *)
112{
113 // Main loop
114 // Called for each event
115
116 Float_t daughter1Mass=AliPID::ParticleMass(AliPID::kKaon);
117 Float_t daughter2Mass=AliPID::ParticleMass(AliPID::kPion);
118
119 if (!fESD) {
120 Printf("ERROR: fESD not available");
121 return;
122 }
123
124 const AliESDVertex* vertex = GetEventVertex(fESD);
125 if(!vertex) return;
126
127 Double_t ptrackpos[3], ptrackneg[3];
128
129 TLorentzVector p4pos;
130 TLorentzVector p4neg;
131 TLorentzVector p4comb;
132
133 for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
134 AliESDtrack* trackpos = fESD->GetTrack(iTracks);
135 if (!trackpos) {
136 Printf("ERROR: Could not receive track %d", iTracks);
137 continue;
138 }
139
140 if (trackpos->GetSign() > 0) continue;
141
142 trackpos->GetPxPyPz(ptrackpos);
143
144 Float_t nSigmaToVertex = GetSigmaToVertex(trackpos);
145
146 Float_t bpos[2];
147 Float_t bCovpos[3];
148 trackpos->GetImpactParameters(bpos,bCovpos);
149
150 if (bCovpos[0]<=0 || bCovpos[2]<=0) {
151 Printf("Estimated b resolution lower or equal zero!");
152 bCovpos[0]=0; bCovpos[2]=0;
153 }
154
155 Float_t dcaToVertexXYpos = bpos[0];
156 Float_t dcaToVertexZpos = bpos[1];
157
158 if(nSigmaToVertex>=4) continue;
159 if((dcaToVertexXYpos>3.0)||(dcaToVertexZpos>3.0)) continue;
160
161 TVector3 posTrackMom(ptrackpos[0],ptrackpos[1],ptrackpos[2]);
162
163 if(posTrackMom.Perp()<=0.25) continue;
164
165 //Uncomment the following block if the Like Sign is made of K- kink + negative track
166
167 //Int_t indexKink=trackpos->GetKinkIndex(0);
168 //Int_t kaonKinkFlag=0;
169// if(indexKink<0){
170//
171// AliESDkink *kink=fESD->GetKink(TMath::Abs(IndexKink)-1);
172// const TVector3 motherMfromKink(kink->GetMotherP());
173// const TVector3 daughterMKink(kink->GetDaughterP());
174// Float_t qt=kink->GetQt();
175//
176// Double_t maxDecAngKmu=f1->Eval(motherMfromKink.Mag(),0.,0.,0.);
177// Double_t maxDecAngpimu=f2->Eval(motherMfromKink.Mag(),0.,0.,0.);
178//
179// Float_t kinkAngle=TMath::RadToDeg()*kink->GetAngle(2);
180//
181// Float_t energyDaughterMu=TMath::Sqrt(daughterMKink.Mag()*daughterMKink.Mag()+0.105658*0.105658);
182// Float_t p1XM= motherMfromKink.Px();
183// Float_t p1YM= motherMfromKink.Py();
184// Float_t p1ZM= motherMfromKink.Pz();
185// Float_t p2XM= daughterMKink.Px();
186// Float_t p2YM= daughterMKink.Py();
187// Float_t p2ZM= daughterMKink.Pz();
188// Float_t p3Daughter=TMath::Sqrt(((p1XM-p2XM)*(p1XM-p2XM))+((p1YM-p2YM)*(p1YM-p2YM))+((p1ZM-p2ZM)*(p1ZM-p2ZM)));
189// Double_t invariantMassKmu= TMath::Sqrt((energyDaughterMu+p3Daughter)*(energyDaughterMu+p3Daughter)-motherMfromKink.Mag()*motherMfromKink.Mag());
190//
191// if((kinkAngle>maxDecAngpimu)&&(qt>0.05)&&(qt<0.25)&&((kink->GetR()>110.)&&(kink->GetR()<230.))&&(TMath::Abs(posTrackMom.Eta())<1.1)&&(invariantMassKmu<0.6)) {
192//
193// if(posTrackMom.Mag()<=1.1) {
194// kaonKinkFlag=1;
195// }
196// else
197// if (kinkAngle<maxDecAngKmu) {
198// kaonKinkFlag=1;
199// }
200// }
201//
202// } //End Kink Information
203//
204// if(kaonKinkFlag==0) continue;
205// if(kaonKinkFlag==1) p4pos.SetVectM(posTrackMom,daughter1Mass);
206
207 // Comment the following statements till the "for" if the Like Sign of K- kink + negative track is needed
208
209 UInt_t status=trackpos->GetStatus();
210 if((status&AliESDtrack::kTPCrefit)==0) continue;
211 if(trackpos->GetTPCclusters(0)<50) continue;
212 if((trackpos->GetTPCchi2()/trackpos->GetTPCclusters(0))>3.5) continue;
213 Double_t extCovPos[15];
214 trackpos->GetExternalCovariance(extCovPos);
215 if(extCovPos[0]>2) continue;
216 if(extCovPos[2]>2) continue;
217 if(extCovPos[5]>0.5) continue;
218 if(extCovPos[9]>0.5) continue;
219 if(extCovPos[14]>2) continue;
220
221 p4pos.SetVectM(posTrackMom,daughter1Mass);
222
223 for (Int_t j=0; j<fESD->GetNumberOfTracks(); j++) {
224 if(iTracks==j) continue;
225 AliESDtrack* trackneg=fESD->GetTrack(j);
226 if (trackneg->GetSign() > 0) continue;
227
228 trackneg->GetPxPyPz(ptrackneg);
229 Float_t negSigmaToVertex = GetSigmaToVertex(trackneg);
230
231 Float_t bneg[2];
232 Float_t bCovneg[3];
233 trackneg->GetImpactParameters(bneg,bCovneg);
234 if (bCovneg[0]<=0 || bCovneg[2]<=0) {
235 Printf("Estimated b resolution lower or equal zero!");
236 bCovneg[0]=0; bCovneg[2]=0;
237 }
238
239 Float_t dcaToVertexXYneg = bneg[0];
240 Float_t dcaToVertexZneg = bneg[1];
241
242 if(negSigmaToVertex>=4) continue;
243 if((dcaToVertexXYneg>3.0)||(dcaToVertexZneg>3.0)) continue;
244
245 TVector3 negTrackMom(ptrackneg[0],ptrackneg[1],ptrackneg[2]);
246
247 if(negTrackMom.Perp()<=0.25) continue;
248
249 UInt_t statusneg=trackneg->GetStatus();
250
251 if((statusneg&AliESDtrack::kTPCrefit)==0) continue;
252
253 if(trackneg->GetTPCclusters(0)<50) continue;
254 if((trackneg->GetTPCchi2()/trackneg->GetTPCclusters(0))>3.5) continue;
255 Double_t extCovneg[15];
256 trackneg->GetExternalCovariance(extCovneg);
257 if(extCovneg[0]>2) continue;
258 if(extCovneg[2]>2) continue;
259 if(extCovneg[5]>0.5) continue;
260 if(extCovneg[9]>0.5) continue;
261 if(extCovneg[14]>2) continue;
262
263 p4neg.SetVectM(negTrackMom, daughter2Mass);
264
265 p4comb=p4pos;
266 p4comb+=p4neg;
267 fNegKaonLikeSign->Fill(p4comb.M());
268
269 } //inner track loop
270
271 } //outer track loop
272
273 PostData(1, fListOfHistos);
274}
275
276//________________________________________________________________________
277void AliResonanceKinkLikeSign::Terminate(Option_t *)
278{
279
280}
281
282//____________________________________________________________________//
283
284Float_t AliResonanceKinkLikeSign::GetSigmaToVertex(AliESDtrack* esdTrack) const {
285 // Calculates the number of sigma to the vertex.
286
287 Float_t b[2];
288 Float_t bRes[2];
289 Float_t bCov[3];
290
291 esdTrack->GetImpactParameters(b,bCov);
292
293 if (bCov[0]<=0 || bCov[2]<=0) {
294 //AliDebug(1, "Estimated b resolution lower or equal zero!");
295 bCov[0]=0; bCov[2]=0;
296 }
297 bRes[0] = TMath::Sqrt(bCov[0]);
298 bRes[1] = TMath::Sqrt(bCov[2]);
299
300 if (bRes[0] == 0 || bRes[1] ==0) return -1;
301
302 Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
303
304 if (TMath::Exp(-d * d / 2) < 1e-10) return 1000;
305
306 d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
307
308 return d;
309}
310
311//____________________________________________________________________//
312
313const AliESDVertex* AliResonanceKinkLikeSign::GetEventVertex(const AliESDEvent* esd) const
314
315{
316 // Get the vertex
317
318 const AliESDVertex* vertex = esd->GetPrimaryVertex();
319
320 if((vertex->GetStatus()==kTRUE)&&(vertex->GetNContributors()>2)) return vertex;
321 else
322 {
323 vertex = esd->GetPrimaryVertexSPD();
324 if((vertex->GetStatus()==kTRUE)&&(vertex->GetNContributors()>2)) return vertex;
325 else
326 return 0;
327 }
328}
329