]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - MUON/MUONmassPlot_ESD.C
- Compute parameter covariances including absorber dispersion effects
[u/mrichter/AliRoot.git] / MUON / MUONmassPlot_ESD.C
index 9f1887a2230b7c8d38707a049e740ec82305f364..71e405e94058c3060ee1b8367b42b9c8d59c7b09 100644 (file)
@@ -1,23 +1,33 @@
 #if !defined(__CINT__) || defined(__MAKECINT__)
 // ROOT includes
+#include "TTree.h"
 #include "TBranch.h"
 #include "TClonesArray.h"
 #include "TLorentzVector.h"
 #include "TFile.h"
 #include "TH1.h"
+#include "TH2.h"
 #include "TParticle.h"
 #include "TTree.h"
 #include <Riostream.h>
+#include <TGeoManager.h>
+#include <TROOT.h>
 
 // STEER includes
 #include "AliRun.h"
+#include "AliLog.h"
 #include "AliRunLoader.h"
 #include "AliHeader.h"
 #include "AliLoader.h"
 #include "AliStack.h"
-#include "AliESD.h"
+#include "AliMagFMaps.h"
+#include "AliESDEvent.h"
+#include "AliESDVertex.h"
+#include "AliTracker.h"
 
 // MUON includes
+#include "AliMUONTrackParam.h"
+#include "AliMUONTrackExtrap.h"
 #include "AliESDMuonTrack.h"
 #endif
 //
 // using Invariant Mass for rapidity.
 
 // Arguments:
+//   ExtrapToVertex (default -1)
+//     <0: no extrapolation;
+//     =0: extrapolation to (0,0,0);
+//     >0: extrapolation to ESDVertex if available, else to (0,0,0)
 //   FirstEvent (default 0)
 //   LastEvent (default 0)
 //   ResType (default 553)
@@ -50,8 +64,8 @@
 
 // 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, 
+Bool_t MUONmassPlot(char* filename = "generated/galice.root", Int_t ExtrapToVertex = -1, char* geoFilename = "geometry.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)
 {
@@ -69,14 +83,19 @@ Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t
   //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 *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", 100, 0., 20.);
+  TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "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 *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.);
+  TH2F *hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",480,0.,12.,80,0.,20.);
+  TH2F *hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",480,0.,12.,80,0.,20.);
   TH1F *hInvMassRes;
+  TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
 
   if (ResType == 553) {
     hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
@@ -88,40 +107,61 @@ Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t
   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.);
+  TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
+  TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
 
 
   // settings
   Int_t EventInMass = 0;
+  Int_t EventInMassMatch = 0;
+  Int_t NbTrigger = 0;
+
   Float_t muonMass = 0.105658389;
 //   Float_t UpsilonMass = 9.46037;
 //   Float_t JPsiMass = 3.097;
 
-  Double_t thetaX, thetaY, pYZ;
+  Int_t fCharge1, fCharge2;
   Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
   Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
-  Int_t fCharge, fCharge2;
 
   Int_t ntrackhits, nevents;
   Double_t fitfmin;
-
+  Double_t fZVertex=0;
+  Double_t fYVertex=0;
+  Double_t fXVertex=0;
+  Double_t errXVtx=0;
+  Double_t errYVtx=0;
  
   TLorentzVector fV1, fV2, fVtot;
+
+  // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
+  if (!gGeoManager) {
+    TGeoManager::Import(geoFilename);
+    if (!gGeoManager) {
+      Error("MUONmass_ESD", "getting geometry from file %s failed", filename);
+      return kFALSE;
+    }
+  }
   
+  // set mag field
+  // waiting for mag field in CDB 
+  printf("Loading field map...\n");
+  AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
+  AliTracker::SetFieldMap(field, kFALSE);
+
   // 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);
+    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);
@@ -130,46 +170,72 @@ Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t
     return kFALSE;
   }
   
+  AliESDEvent* esd = new AliESDEvent();
+  TTree* tree = (TTree*) esdFile->Get("esdTree");
+  if (!tree) {
+    Error("CheckESD", "no ESD tree found");
+    return kFALSE;
+  }
+//  tree->SetBranchAddress("ESD", &esd);
+  esd->ReadFromTree(tree);
+  
+
   runLoader->LoadHeader();
   nevents = runLoader->GetNumberOfEvents();
-        
+  
+  AliMUONTrackParam trackParam;
+
   // 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);
+    tree->GetEvent(iEvent);
     if (!esd) {
-      Error("MUONmass_ESD", "no ESD object found for event %d", iEvent);
+      Error("CheckESD", "no ESD object found for event %d", iEvent);
       return kFALSE;
     }
 
-    Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; //
+    // get the SPD reconstructed vertex (vertexer) and fill the histogram
+    AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
+    if (Vertex->GetNContributors()) {
+      fZVertex = Vertex->GetZv();
+      fYVertex = Vertex->GetYv();
+      fXVertex = Vertex->GetXv();
+      errXVtx = Vertex->GetXRes();
+      errYVtx = Vertex->GetYRes();
+    }
+    hPrimaryVertex->Fill(fZVertex);
+
+    Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
 
     //    printf("\n Nb of events analysed: %d\r",iEvent);
-    //   cout << " number of tracks: " << nrectracks  <<endl;
+    //      cout << " number of tracks: " << nTracks  <<endl;
   
+    // set the magnetic field for track extrapolations
+    AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
     // 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);
-
+      AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack)));
+
+      // extrapolate to vertex if required and available
+      if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
+        trackParam.GetParamFromUncorrected(*muonTrack);
+       AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
+       trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
+      } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
+        trackParam.GetParamFromUncorrected(*muonTrack);
+       AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
+       trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
+      }
+
+      fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
+      
+      muonTrack->LorentzP(fV1);
+      
       ntrackhits = muonTrack->GetNHit();
       fitfmin    = muonTrack->GetChi2();
 
@@ -185,7 +251,7 @@ Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t
       // 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);
+//          fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge1);
 
       // condition for good track (Chi2Cut and PtCut)
 
@@ -196,26 +262,35 @@ Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t
        hPMuon->Fill(p1);
        hChi2PerDof->Fill(ch1);
        hRapMuon->Fill(rapMuon1);
-
+       if (fCharge1 > 0) {
+         hPtMuonPlus->Fill(pt1);
+         hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
+       } else {
+         hPtMuonMinus->Fill(pt1);
+         hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
+       }
        // 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()));
+         AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2)));
+         
+         // extrapolate to vertex if required and available
+         if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
+           trackParam.GetParamFromUncorrected(*muonTrack2);
+           AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
+           trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
+         } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
+            trackParam.GetParamFromUncorrected(*muonTrack2);
+           AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
+           trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
+         }
+         
+         fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));
 
-         fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
-         fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
+         muonTrack2->LorentzP(fV2);
 
-         ntrackhits = muonTrack->GetNHit();
-         fitfmin    = muonTrack->GetChi2();
+         ntrackhits = muonTrack2->GetNHit();
+         fitfmin    = muonTrack2->GetChi2();
 
          // transverse momentum
          Float_t pt2 = fV2.Pt();
@@ -227,7 +302,7 @@ Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t
          if ((ch2 < Chi2Cut) && (pt2 > PtCutMin)  && (pt2 < PtCutMax)) {
 
            // condition for opposite charges
-           if ((fCharge * fCharge2) == -1) {
+           if ((fCharge1 * fCharge2) == -1) {
 
              // invariant mass
              fVtot = fV1 + fV2;
@@ -236,35 +311,77 @@ Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t
              // fill histos hInvMassAll and hInvMassRes
              hInvMassAll->Fill(invMass);
              hInvMassRes->Fill(invMass);
-
+             hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
+             Int_t ptTrig;
+             if (ResType == 553) 
+               ptTrig =  0x20;// mask for Hpt unlike sign pair
+             else 
+               ptTrig =  0x10;// mask for Lpt unlike sign pair
+
+             if (esd->GetTriggerMask() &  ptTrig) NbTrigger++; 
              if (invMass > massMin && invMass < massMax) {
                EventInMass++;
+               if (muonTrack2->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig))// match with trigger
+                 EventInMassMatch++;
+
                hRapResonance->Fill(fVtot.Rapidity());
                hPtResonance->Fill(fVtot.Pt());
              }
 
-           } // if (fCharge * fCharge2) == -1)
+           } // if (fCharge1 * fCharge2) == -1)
          } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
+         delete muonTrack2;
        } //  for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
       } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
+      delete muonTrack;
     } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
 
     hNumberOfTrack->Fill(nTracks);
+    //    esdFile->Delete();
   } // for (Int_t iEvent = FirstEvent;
 
+// Loop over events for bg event
+
+  Double_t thetaPlus,  phiPlus;
+  Double_t thetaMinus, phiMinus;
+  Float_t PtMinus, PtPlus;
+  
+  for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {
+
+    hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
+    hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
+    PtPlus = hPtMuonPlus->GetRandom();
+    PtMinus = hPtMuonMinus->GetRandom();
+
+    fPxRec1  = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
+    fPyRec1  = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
+    fPzRec1  = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
+
+    fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
+    fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
+
+    fPxRec2  = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
+    fPyRec2  = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
+    fPzRec2  = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
+
+    fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
+    fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
+
+    // invariant mass
+    fVtot = fV1 + fV2;
+      
+    // fill histos hInvMassAll and hInvMassRes
+    hInvMassBg->Fill(fVtot.M());
+    hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
+  }
+
   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 << endl;
   cout << "EventInMass " << EventInMass << endl;
+  cout << "NbTrigger " << NbTrigger << endl;
+  cout << "EventInMass match with trigger " << EventInMassMatch << endl;
 
   return kTRUE;
 }