]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONmassPlot.C
Add ResetDecayTable() and SsetDecayTable() methods.
[u/mrichter/AliRoot.git] / MUON / MUONmassPlot.C
1 // macro to make invariant mass plots
2 // for combinations of 2 muons with opposite charges,
3 // from root file "MUONtrackReco.root" containing the result of track reconstruction,
4 // generated by the macro "MUONrecoNtuple.C".
5 // Histograms are stored on the "MUONmassPlot.root" file.
6 // A model for macros using the Ntuple in the file "MUONtrackReco.root"
7 // may be found in "MUONtrackRecoModel.C":
8 // it has been obtained by reading with Root a file "MUONtrackReco.root",
9 // and executing the command:
10 // MUONtrackReco->MakeCode("MUONtrackRecoModel.C")
11
12 // Arguments:
13 //   FirstEvent (default 0)
14 //   LastEvent (default 0)
15 //   ResType (default 553)
16 //      553 for Upsilon, anything else for J/Psi
17 //   NSigma (default 3)
18 //      the number of combinations is counted around the resonance mass
19 //      within +/- NSigma times the nominal sigma's
20 //      (0.099 GeV for Upsilon, 0.0615 GeV for J/Psi)
21 //   Chi2Cut (default 100)
22 //      to keep only tracks with chi2 per d.o.f. < Chi2Cut
23 //   PtCut (default 1)
24 //      to keep only tracks with transverse momentum > PtCut
25
26 void MUONmassPlot(Int_t FirstEvent = 0, Int_t LastEvent = 0, Int_t ResType = 553, Float_t Nsig = 3., Float_t Chi2Cut = 100., Float_t PtCut = 1.)
27 {
28   cout << "MUONmassPlot" << endl;
29   cout << "FirstEvent" << FirstEvent << endl;
30   cout << "LastEvent" << LastEvent << endl;
31   cout << "ResType" << ResType << endl;
32   cout << "Nsig" << Nsig << endl;
33   cout << "Chi2Cut" << Chi2Cut << endl;
34   cout << "PtCut" << PtCut << endl;
35
36   //////////////////////////////////////////////////////////
37   //   This file has been automatically generated 
38   //     (Thu Sep 21 14:53:11 2000 by ROOT version2.25/02)
39   //   from TTree MUONtrackReco/MUONtrackReco
40   //   found on file: MUONtrackReco.root
41   //////////////////////////////////////////////////////////
42
43
44   //Reset ROOT and connect tree file
45   gROOT->Reset();
46   TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject("MUONtrackReco.root");
47   if (!f) {
48     f = new TFile("MUONtrackReco.root");
49   }
50   TTree *MUONtrackReco = (TTree*)gDirectory->Get("MUONtrackReco");
51
52 //Declaration of leaves types
53   Int_t           fEvent;
54   UInt_t          fUniqueID;
55   UInt_t          fBits;
56   Int_t           Tracks_;
57   Int_t           Tracks_fCharge[5];
58   Float_t         Tracks_fPxRec[5];
59   Float_t         Tracks_fPyRec[5];
60   Float_t         Tracks_fPzRec[5];
61   Float_t         Tracks_fZRec[5];
62   Float_t         Tracks_fZRec1[5];
63   Int_t           Tracks_fNHits[5];
64   Float_t         Tracks_fChi2[5];
65   Float_t         Tracks_fPxGen[5];
66   Float_t         Tracks_fPyGen[5];
67   Float_t         Tracks_fPzGen[5];
68   UInt_t          Tracks_fUniqueID[5];
69   UInt_t          Tracks_fBits[5];
70
71   //Set branch addresses
72   //MUONtrackReco->SetBranchAddress("Header",&Header);
73   MUONtrackReco->SetBranchAddress("fEvent",&fEvent);
74   MUONtrackReco->SetBranchAddress("fUniqueID",&fUniqueID);
75   MUONtrackReco->SetBranchAddress("fBits",&fBits);
76   MUONtrackReco->SetBranchAddress("Tracks_",&Tracks_);
77   MUONtrackReco->SetBranchAddress("Tracks.fCharge",Tracks_fCharge);
78   MUONtrackReco->SetBranchAddress("Tracks.fPxRec",Tracks_fPxRec);
79   MUONtrackReco->SetBranchAddress("Tracks.fPyRec",Tracks_fPyRec);
80   MUONtrackReco->SetBranchAddress("Tracks.fPzRec",Tracks_fPzRec);
81   MUONtrackReco->SetBranchAddress("Tracks.fZRec",Tracks_fZRec);
82   MUONtrackReco->SetBranchAddress("Tracks.fZRec1",Tracks_fZRec1);
83   MUONtrackReco->SetBranchAddress("Tracks.fNHits",Tracks_fNHits);
84   MUONtrackReco->SetBranchAddress("Tracks.fChi2",Tracks_fChi2);
85   MUONtrackReco->SetBranchAddress("Tracks.fPxGen",Tracks_fPxGen);
86   MUONtrackReco->SetBranchAddress("Tracks.fPyGen",Tracks_fPyGen);
87   MUONtrackReco->SetBranchAddress("Tracks.fPzGen",Tracks_fPzGen);
88   MUONtrackReco->SetBranchAddress("Tracks.fUniqueID",Tracks_fUniqueID);
89   MUONtrackReco->SetBranchAddress("Tracks.fBits",Tracks_fBits);
90
91 //     This is the loop skeleton
92 //       To read only selected branches, Insert statements like:
93 // MUONtrackReco->SetBranchStatus("*",0);  // disable all branches
94 // TTreePlayer->SetBranchStatus("branchname",1);  // activate branchname
95
96   Int_t nentries = MUONtrackReco->GetEntries();
97
98   Int_t nbytes = 0;
99   //   for (Int_t i=0; i<nentries;i++) {
100   //      nbytes += MUONtrackReco->GetEntry(i);
101   //   }
102
103   /////////////////////////////////////////////////////////////////
104   // Here comes the specialized part for MUONmassPlot
105   /////////////////////////////////////////////////////////////////
106
107   // File for histograms and histogram booking
108   TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
109   TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
110   TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
111   TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 240, 0., 12.);
112   if (ResType = 553) TH1F *hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
113   else TH1F *hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 1., 5.);
114
115   // Loop over events
116   for (Int_t event = FirstEvent; event <= TMath::Min(LastEvent, nentries - 1); event++) {
117     // get current event
118     nbytes += MUONtrackReco->GetEntry(event);
119     // loop over all reconstructed tracks (also first track of combination)
120     for (Int_t t1 = 0; t1 < Tracks_; t1++) {
121       // transverse momentum
122       Float_t pt1 = TMath::Sqrt(Tracks_fPxRec[t1] * Tracks_fPxRec[t1] + Tracks_fPyRec[t1] * Tracks_fPyRec[t1]);
123       // chi2 per d.o.f.
124       Float_t ch1 = Tracks_fChi2[t1] / (2.0 * Tracks_fNHits[t1] - 5);
125       // condition for good track (Chi2Cut and PtCut)
126       if ((ch1 < Chi2Cut) && (pt1 > PtCut)) {
127         // fill histos hPtMuon and hChi2PerDof
128         hPtMuon->Fill(pt1);
129         hChi2PerDof->Fill(ch1);
130         // loop over second track of combination
131         for (Int_t t2 = t1 + 1; t2 < Tracks_; t2++) {
132           // transverse momentum
133           Float_t pt2 = TMath::Sqrt(Tracks_fPxRec[t2] * Tracks_fPxRec[t2] + Tracks_fPyRec[t2] * Tracks_fPyRec[t2]);
134           // chi2 per d.o.f.
135           Float_t ch2 = Tracks_fChi2[t2] / (2.0 * Tracks_fNHits[t2] - 5);
136           // condition for good track (Chi2Cut and PtCut)
137           if ((ch2 < Chi2Cut) && (pt2 > PtCut)) {
138             // condition for opposite charges
139             if ((Tracks_fCharge[t1] * Tracks_fCharge[t2]) == -1) {
140               // invariant mass
141               Float_t invMass = MuPlusMuMinusMass(Tracks_fPxRec[t1], Tracks_fPyRec[t1], Tracks_fPzRec[t1], Tracks_fPxRec[t2], Tracks_fPyRec[t2], Tracks_fPzRec[t2]);
142               // fill histos hInvMassAll and hInvMassRes
143               hInvMassAll->Fill(invMass);
144               hInvMassRes->Fill(invMass);
145             } //if ((Tracks_fCharge[t1] * Tracks_fCharge[t2]) == -1)
146           } // if ((Tracks_fChi2[t2] < Chi2Cut) && (pt2 > PtCut))
147         } // for (Int_t t2 = t1 + 1; t2 < Tracks_; t2++)
148       } // if ((Tracks_fChi2[t1] < Chi2Cut) && (pt1 > PtCut))
149     } // for (Int_t t1 = 0; t1 < Tracks_; t1++)
150   } // for (Int_t event = FirstEvent;
151
152   histoFile->Write();
153   histoFile->Close();
154 }
155
156 Float_t MuPlusMuMinusMass(Float_t Px1, Float_t Py1, Float_t Pz1, Float_t Px2, Float_t Py2, Float_t Pz2)
157 {
158   Float_t muonMass = 0.10566;
159   Float_t e1 = TMath::Sqrt(muonMass * muonMass + Px1 * Px1 + Py1 * Py1 + Pz1 * Pz1);
160   Float_t e2 = TMath::Sqrt(muonMass * muonMass + Px2 * Px2 + Py2 * Py2 + Pz2 * Pz2);
161   return (TMath::Sqrt(2.0 * (muonMass * muonMass + e1 * e2 - Px1 * Px2 - Py1 * Py2 - Pz1 * Pz2)));
162 }