New macro for invariant mass analysis from ESD (Christian)
authormartinez <martinez@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 7 Apr 2004 14:56:33 +0000 (14:56 +0000)
committermartinez <martinez@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 7 Apr 2004 14:56:33 +0000 (14:56 +0000)
MUON/MUONmassPlot_ESD.C [new file with mode: 0644]
MUON/README

diff --git a/MUON/MUONmassPlot_ESD.C b/MUON/MUONmassPlot_ESD.C
new file mode 100644 (file)
index 0000000..9f1887a
--- /dev/null
@@ -0,0 +1,271 @@
+#if !defined(__CINT__) || defined(__MAKECINT__)
+// ROOT includes
+#include "TBranch.h"
+#include "TClonesArray.h"
+#include "TLorentzVector.h"
+#include "TFile.h"
+#include "TH1.h"
+#include "TParticle.h"
+#include "TTree.h"
+#include <Riostream.h>
+
+// STEER includes
+#include "AliRun.h"
+#include "AliRunLoader.h"
+#include "AliHeader.h"
+#include "AliLoader.h"
+#include "AliStack.h"
+#include "AliESD.h"
+
+// MUON includes
+#include "AliESDMuonTrack.h"
+#endif
+//
+// Macro MUONmassPlot.C for ESD
+// Ch. Finck, Subatech, April. 2004
+//
+
+// macro to make invariant mass plots
+// for combinations of 2 muons with opposite charges,
+// from root file "MUON.tracks.root" containing the result of track reconstruction.
+// Histograms are stored on the "MUONmassPlot.root" file.
+// introducing TLorentzVector for parameter calculations (Pt, P,rap,etc...)
+// using Invariant Mass for rapidity.
+
+// Arguments:
+//   FirstEvent (default 0)
+//   LastEvent (default 0)
+//   ResType (default 553)
+//      553 for Upsilon, anything else for J/Psi
+//   Chi2Cut (default 100)
+//      to keep only tracks with chi2 per d.o.f. < Chi2Cut
+//   PtCutMin (default 1)
+//      to keep only tracks with transverse momentum > PtCutMin
+//   PtCutMax (default 10000)
+//      to keep only tracks with transverse momentum < PtCutMax
+//   massMin (default 9.17 for Upsilon) 
+//      &  massMax (default 9.77 for Upsilon) 
+//         to calculate the reconstruction efficiency for resonances with invariant mass
+//         massMin < mass < massMax.
+
+// Add parameters and histograms for analysis 
+
+Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t LastEvent = 10000,
+                 char* esdFileName = "AliESDs.root", 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)
+{
+  cout << "MUONmassPlot " << endl;
+  cout << "FirstEvent " << FirstEvent << endl;
+  cout << "LastEvent " << LastEvent << endl;
+  cout << "ResType " << ResType << endl;
+  cout << "Chi2Cut " << Chi2Cut << endl;
+  cout << "PtCutMin " << PtCutMin << endl;
+  cout << "PtCutMax " << PtCutMax << endl;
+  cout << "massMin " << massMin << endl;
+  cout << "massMax " << massMax << endl;
+
+  //Reset ROOT and connect tree file
+  gROOT->Reset();
+
+
+  // File for histograms and histogram booking
+  TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
+  TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
+  TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", 100, 0., 200.);
+  TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
+  TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
+  TH1F *hInvMassRes;
+
+  if (ResType == 553) {
+    hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
+  } else {
+    hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
+  }
+
+  TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
+  TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
+  TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
+  TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
+
+
+  // settings
+  Int_t EventInMass = 0;
+  Float_t muonMass = 0.105658389;
+//   Float_t UpsilonMass = 9.46037;
+//   Float_t JPsiMass = 3.097;
+
+  Double_t thetaX, thetaY, pYZ;
+  Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
+  Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
+  Int_t fCharge, fCharge2;
+
+  Int_t ntrackhits, nevents;
+  Double_t fitfmin;
+
+  TLorentzVector fV1, fV2, fVtot;
+  
+  // open run loader and load gAlice, kinematics and header
+  AliRunLoader* runLoader = AliRunLoader::Open(filename);
+  if (!runLoader) {
+    Error("MUONmass_ESD", "getting run loader from file %s failed", 
+           filename);
+    return kFALSE;
+  }
+
+  runLoader->LoadgAlice();
+  gAlice = runLoader->GetAliRun();
+  if (!gAlice) {
+    Error("MUONmass_ESD", "no galice object found");
+    return kFALSE;
+  }
+  
+
+  // open the ESD file
+  TFile* esdFile = TFile::Open(esdFileName);
+  if (!esdFile || !esdFile->IsOpen()) {
+    Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName);
+    return kFALSE;
+  }
+  
+  runLoader->LoadHeader();
+  nevents = runLoader->GetNumberOfEvents();
+        
+  // Loop over events
+  for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) {
+
+    // get current event
+    runLoader->GetEvent(iEvent);
+    
+    // get the event summary data
+    char esdName[256]; 
+    sprintf(esdName, "ESD%d", iEvent);
+    AliESD* esd = (AliESD*) esdFile->Get(esdName);
+    if (!esd) {
+      Error("MUONmass_ESD", "no ESD object found for event %d", iEvent);
+      return kFALSE;
+    }
+
+    Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; //
+
+    //    printf("\n Nb of events analysed: %d\r",iEvent);
+    //   cout << " number of tracks: " << nrectracks  <<endl;
+  
+    // loop over all reconstructed tracks (also first track of combination)
+    for (Int_t iTrack = 0; iTrack <  nTracks;  iTrack++) {
+
+      AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
+
+      thetaX = muonTrack->GetThetaX();
+      thetaY = muonTrack->GetThetaY();
+
+      pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
+      fPzRec1  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
+      fPxRec1  = fPzRec1 * TMath::Tan(thetaX);
+      fPyRec1  = fPzRec1 * TMath::Tan(thetaY);
+      fCharge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
+
+      fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
+      fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
+
+      ntrackhits = muonTrack->GetNHit();
+      fitfmin    = muonTrack->GetChi2();
+
+      // transverse momentum
+      Float_t pt1 = fV1.Pt();
+
+      // total momentum
+      Float_t p1 = fV1.P();
+
+      // Rapidity
+      Float_t rapMuon1 = fV1.Rapidity();
+
+      // chi2 per d.o.f.
+      Float_t ch1 =  fitfmin / (2.0 * ntrackhits - 5);
+//      printf(" px %f py %f pz %f NHits %d  Norm.chi2 %f charge %d\n", 
+//          fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge);
+
+      // condition for good track (Chi2Cut and PtCut)
+
+      if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
+
+       // fill histos hPtMuon and hChi2PerDof
+       hPtMuon->Fill(pt1);
+       hPMuon->Fill(p1);
+       hChi2PerDof->Fill(ch1);
+       hRapMuon->Fill(rapMuon1);
+
+       // loop over second track of combination
+       for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
+         
+         AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack2);
+
+         thetaX = muonTrack->GetThetaX();
+         thetaY = muonTrack->GetThetaY();
+
+         pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
+         fPzRec2  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
+         fPxRec2  = fPzRec2 * TMath::Tan(thetaX);
+         fPyRec2  = fPzRec2 * TMath::Tan(thetaY);
+         fCharge2 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
+
+         fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
+         fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
+
+         ntrackhits = muonTrack->GetNHit();
+         fitfmin    = muonTrack->GetChi2();
+
+         // transverse momentum
+         Float_t pt2 = fV2.Pt();
+
+         // chi2 per d.o.f.
+         Float_t ch2 = fitfmin  / (2.0 * ntrackhits - 5);
+
+         // condition for good track (Chi2Cut and PtCut)
+         if ((ch2 < Chi2Cut) && (pt2 > PtCutMin)  && (pt2 < PtCutMax)) {
+
+           // condition for opposite charges
+           if ((fCharge * fCharge2) == -1) {
+
+             // invariant mass
+             fVtot = fV1 + fV2;
+             Float_t invMass = fVtot.M();
+                   
+             // fill histos hInvMassAll and hInvMassRes
+             hInvMassAll->Fill(invMass);
+             hInvMassRes->Fill(invMass);
+
+             if (invMass > massMin && invMass < massMax) {
+               EventInMass++;
+               hRapResonance->Fill(fVtot.Rapidity());
+               hPtResonance->Fill(fVtot.Pt());
+             }
+
+           } // if (fCharge * fCharge2) == -1)
+         } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
+       } //  for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
+      } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
+    } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
+
+    hNumberOfTrack->Fill(nTracks);
+  } // for (Int_t iEvent = FirstEvent;
+
+  histoFile->Write();
+  histoFile->Close();
+
+  cout << "MUONmassPlot " << endl;
+  cout << "FirstEvent " << FirstEvent << endl;
+  cout << "LastEvent " << LastEvent << endl;
+  cout << "ResType " << ResType << endl;
+  cout << "Chi2Cut " << Chi2Cut << endl;
+  cout << "PtCutMin " << PtCutMin << endl;
+  cout << "PtCutMax " << PtCutMax << endl;
+  cout << "massMin " << massMin << endl;
+  cout << "massMax " << massMax << endl;
+  cout << "EventInMass " << EventInMass << endl;
+
+  return kTRUE;
+}
+
index adbbbec..85983ac 100644 (file)
@@ -255,42 +255,27 @@ gener->Init();
 ===========================================================
  csh Script for the full reconstruction
 ===========================================================
-aliroot -b << EOF  
-gAlice->Run(NumberOfEvents,"YourConfig.C");
-.q
-EOF
 
-aliroot -b << EOF
-AliRunDigitizer   * manager = new AliRunDigitizer(1,1);
-manager->SetInputStream(0,"galice.root");
-AliMUONDigitizer * dMUON   = new AliMUONSDigitizerv1(manager);
-manager->AddDigitizer(dMUON);
-manager->Exec("deb");
-.q.
-EOF
 
-aliroot -b << EOF
-AliRunDigitizer   * manager = new AliRunDigitizer(1,1);
-manager->SetInputStream(0,"galice.root");
-AliMUONDigitizer* dMUON   = new AliMUONDigitizerv2(manager);
-manager->AddDigitizer(dMUON);
-manager->Exec("deb");
+aliroot -b << EOF  
+AliSimulation MuonSim("YourConfig.C")
+MuonSim.Run(XXX)
 .q
 EOF
 
 aliroot -b << EOF
-.includepath $ALICE_ROOT/STEER
-.includepath $ALICE_ROOT/MUON
-.L $ALICE_ROOT/MUON/MUONChallengeTest.C++
-MUONRecoTest("galice.root");
+AliReconstruction MuonRec("galice.root")
+MuonRec.SetRunVertexFinder(kFALSE)
+MuonRec.SetRunLocalReconstruction("MUON")
+MuonRec.SetFillESD("MUON")
+MuonRec.Run()
 .q
 EOF
 
-
 aliroot -b << EOF  
 .includepath $ALICE_ROOT/STEER
 .includepath $ALICE_ROOT/MUON
-.L $ALICE_ROOT/MUON/MUONmassPlot_NewIO.C++
+.L $ALICE_ROOT/MUON/MUONmassPlot_ESD.C++
 MUONmassPlot("galice.root",0,99999);
 .q
 EOF