1 #if !defined(__CINT__) || defined(__MAKECINT__) 5 #include "TClonesArray.h" 6 #include "TLorentzVector.h" 10 #include "TParticle.h" 12 #include <Riostream.h> 13 #include <TGeoManager.h> 47 Bool_t
MUONmassPlot(
const char* esdFileName =
"AliESDs.root",
const char* geoFilename =
"geometry.root",
48 const char* ocdbPath =
"local://$ALICE_ROOT/OCDB",
49 Int_t FirstEvent = 0, Int_t LastEvent = -1, Int_t ExtrapToVertex = -1,
50 Int_t ResType = 553, Float_t Chi2Cut = 100., Float_t PtCutMin = 1.,
51 Float_t PtCutMax = 10000., Float_t massMin = 9.17,Float_t massMax = 9.77)
68 cout <<
"MUONmassPlot " << endl;
69 cout <<
"FirstEvent " << FirstEvent << endl;
70 cout <<
"LastEvent " << ((LastEvent>=0) ? Form(
"%d",LastEvent) :
"all") << endl;
71 cout <<
"ResType " << ResType << endl;
72 cout <<
"Chi2Cut " << Chi2Cut << endl;
73 cout <<
"PtCutMin " << PtCutMin << endl;
74 cout <<
"PtCutMax " << PtCutMax << endl;
75 cout <<
"massMin " << massMin << endl;
76 cout <<
"massMax " << massMax << endl;
83 TFile *histoFile =
new TFile(
"MUONmassPlot.root",
"RECREATE");
84 TH1F *hPtMuon =
new TH1F(
"hPtMuon",
"Muon Pt (GeV/c)", 100, 0., 20.);
85 TH1F *hPtMuonPlus =
new TH1F(
"hPtMuonPlus",
"Muon+ Pt (GeV/c)", 100, 0., 20.);
86 TH1F *hPtMuonMinus =
new TH1F(
"hPtMuonMinus",
"Muon- Pt (GeV/c)", 100, 0., 20.);
87 TH1F *hPMuon =
new TH1F(
"hPMuon",
"Muon P (GeV/c)", 100, 0., 200.);
88 TH1F *hChi2PerDof =
new TH1F(
"hChi2PerDof",
"Muon track chi2/d.o.f.", 100, 0., 20.);
89 TH1F *hInvMassAll =
new TH1F(
"hInvMassAll",
"Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
90 TH1F *hInvMassBg =
new TH1F(
"hInvMassBg",
"Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.);
91 TH2F *hInvMassAll_vs_Pt =
new TH2F(
"hInvMassAll_vs_Pt",
"hInvMassAll_vs_Pt",480,0.,12.,80,0.,20.);
92 TH2F *hInvMassBgk_vs_Pt =
new TH2F(
"hInvMassBgk_vs_Pt",
"hInvMassBgk_vs_Pt",480,0.,12.,80,0.,20.);
94 TH1F *hPrimaryVertex =
new TH1F(
"hPrimaryVertex",
"SPD reconstructed Z vertex",150,-15,15);
97 hInvMassRes =
new TH1F(
"hInvMassRes",
"Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
99 hInvMassRes =
new TH1F(
"hInvMassRes",
"Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
102 TH1F *hNumberOfTrack =
new TH1F(
"hNumberOfTrack",
"nb of track /evt ",20,-0.5,19.5);
103 TH1F *hRapMuon =
new TH1F(
"hRapMuon",
" Muon Rapidity",50,-4.5,-2);
104 TH1F *hRapResonance =
new TH1F(
"hRapResonance",
" Resonance Rapidity",50,-4.5,-2);
105 TH1F *hPtResonance =
new TH1F(
"hPtResonance",
"Resonance Pt (GeV/c)", 100, 0., 20.);
106 TH2F *hThetaPhiPlus =
new TH2F(
"hThetaPhiPlus",
"Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
107 TH2F *hThetaPhiMinus =
new TH2F(
"hThetaPhiMinus",
"Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
111 Int_t EventInMass = 0;
112 Int_t EventInMassMatch = 0;
115 Float_t muonMass = 0.105658389;
119 Int_t fCharge1, fCharge2;
120 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
121 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
131 TLorentzVector fV1, fV2, fVtot;
135 TGeoManager::Import(geoFilename);
137 Error(
"MUONmass_ESD",
"getting geometry from file %s failed", geoFilename);
144 if (!esdFile || !esdFile->IsOpen()) {
145 Error(
"MUONmass_ESD",
"opening ESD file %s failed", esdFileName);
149 TTree*
tree = (TTree*) esdFile->Get(
"esdTree");
151 Error(
"MUONmass_ESD",
"no ESD tree found");
157 if (tree->GetEvent(0) <= 0) {
158 Error(
"MUONmass_ESD",
"no ESD object found for event 0");
176 Int_t nevents = (LastEvent >= 0) ? TMath::Min(LastEvent, (Int_t)tree->GetEntries()-1) : (Int_t)tree->GetEntries()-1;
177 for (Int_t iEvent = FirstEvent; iEvent <= nevents; iEvent++) {
180 if (tree->GetEvent(iEvent) <= 0) {
181 Error(
"MUONmass_ESD",
"no ESD object found for event %d", iEvent);
188 fZVertex = Vertex->
GetZ();
189 fYVertex = Vertex->
GetY();
190 fXVertex = Vertex->
GetX();
194 hPrimaryVertex->Fill(fZVertex);
202 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
214 }
else if ((ExtrapToVertex > 0 && !Vertex->
GetNContributors()) || ExtrapToVertex == 0){
225 Double_t dAbs1 = TMath::Sqrt(xAbs*xAbs + yAbs*yAbs);
232 ntrackhits = muonTrack->
GetNHit();
233 fitfmin = muonTrack->
GetChi2();
236 Float_t pt1 = fV1.Pt();
239 Float_t p1 = fV1.P();
242 Float_t rapMuon1 = fV1.Rapidity();
245 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
256 hChi2PerDof->Fill(ch1);
257 hRapMuon->Fill(rapMuon1);
259 hPtMuonPlus->Fill(pt1);
260 hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
262 hPtMuonMinus->Fill(pt1);
263 hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
266 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
278 }
else if ((ExtrapToVertex > 0 && !Vertex->
GetNContributors()) || ExtrapToVertex == 0){
289 Double_t dAbs2 = TMath::Sqrt(xAbs*xAbs + yAbs*yAbs);
296 ntrackhits = muonTrack2->
GetNHit();
297 fitfmin = muonTrack2->
GetChi2();
300 Float_t pt2 = fV2.Pt();
303 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
318 if ((fCharge1 * fCharge2) == -1) {
322 Float_t invMass = fVtot.M();
325 hInvMassAll->Fill(invMass);
326 hInvMassRes->Fill(invMass);
327 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
335 if (invMass > massMin && invMass < massMax) {
340 hRapResonance->Fill(fVtot.Rapidity());
341 hPtResonance->Fill(fVtot.Pt());
352 hNumberOfTrack->Fill(nTracks);
358 Double_t thetaPlus, phiPlus;
359 Double_t thetaMinus, phiMinus;
360 Float_t PtMinus, PtPlus;
362 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {
364 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
365 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
366 PtPlus = hPtMuonPlus->GetRandom();
367 PtMinus = hPtMuonMinus->GetRandom();
369 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
370 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
371 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
373 fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
374 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
376 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
377 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
378 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
380 fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
381 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
387 hInvMassBg->Fill(fVtot.M());
388 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
395 cout <<
"EventInMass " << EventInMass << endl;
396 cout <<
"NbTrigger " << NbTrigger << endl;
397 cout <<
"EventInMass match with trigger " << EventInMassMatch << endl;
static Double_t AbsZEnd()
Return z-position of absorber end.
TFile * Open(const char *filename, Long64_t &nevents)
AliESDMuonTrack * GetMuonTrack(Int_t i)
UChar_t GetNHit(void) const
Class to describe the MUON tracks in the Event Summary Data class.
virtual Double_t GetY() const
Double_t GetBendingCoor() const
return bending coordinate (cm)
static void SetParamAtVertex(const AliMUONTrackParam &trackParam, AliESDMuonTrack &esdTrack)
Track parameters in ALICE dimuon spectrometer.
Bool_t ContainTrackerData() const
static void GetParamAtFirstCluster(const AliESDMuonTrack &esdTrack, AliMUONTrackParam &trackParam)
ULong64_t GetTriggerMask() const
void SetSpecificStorage(const char *calibType, const char *dbString, Int_t version=-1, Int_t subVersion=-1)
Int_t GetNumberOfMuonTracks() const
Double_t GetInverseBendingMomentum(void) const
Int_t GetRunNumber() const
TGeoManager * gGeoManager
virtual Double_t GetZ() const
Double_t GetChi2(void) const
void LorentzP(TLorentzVector &vP) const
virtual Double_t GetX() const
Int_t GetMatchTrigger() const
void SetDefaultStorage(const char *dbString)
Bool_t MUONmassPlot(const char *esdFileName="AliESDs.root", const char *geoFilename="geometry.root", const char *ocdbPath="local://$ALICE_ROOT/OCDB", Int_t FirstEvent=0, Int_t LastEvent=-1, Int_t ExtrapToVertex=-1, Int_t ResType=553, Float_t Chi2Cut=100., Float_t PtCutMin=1., Float_t PtCutMax=10000., Float_t massMin=9.17, Float_t massMax=9.77)
void ReadFromTree(TTree *tree, Option_t *opt="")
const AliESDVertex * GetVertex() const
virtual Int_t GetNContributors() const
Double_t GetNonBendingCoor() const
return non bending coordinate (cm)
static AliCDBManager * Instance(TMap *entryCache=NULL, Int_t run=-1)