]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
New MFT Analysis Tools
authorauras <antonio.uras@cern.ch>
Sat, 25 Oct 2014 17:57:54 +0000 (19:57 +0200)
committerauras <antonio.uras@cern.ch>
Sat, 25 Oct 2014 17:57:54 +0000 (19:57 +0200)
MFT/AliAnalysisTaskMFTExample.cxx
MFT/AliMFTAnalysisTools.cxx
MFT/AliMFTAnalysisTools.h
MFT/AliMFTPlane.cxx

index f6b10791016ebac983caeaad35a83ee8cb7eaf9a..0a6c54f03c5eb8e8c3279cec87c9640f1d64a344 100755 (executable)
@@ -197,30 +197,38 @@ void AliAnalysisTaskMFTExample::UserExec(Option_t *) {
       recMuon2 = aodEv->GetTrack(jTrack);
       
       AliAODDimuon *dimuon = new AliAODDimuon(recMuon1, recMuon2);
-      if (dimuon->Charge()) continue;
+      if (dimuon->Charge()) {
+       delete dimuon;
+       continue;
+      }
 
       // pt and mass all OS dimuons
       fHistPtDimuonsOS   -> Fill(dimuon->Pt());
       fHistMassDimuonsOS -> Fill(dimuon->Mass());
 
+      delete dimuon;
+
       // pt and mass J/psi dimuons
       if (!isMuon1FromJpsi) continue;
       if (recMuon2->GetLabel() >= 0) {
        mcMuon2 = (AliAODMCParticle*) stackMC->At(recMuon2->GetLabel());
        if (mcMuon2) {
          if (mcMuon2->GetMother() == mcMuon1->GetMother()) {
-           AliAODDimuon *dimuon = new AliAODDimuon;
-           dimuon->SetMuons(recMuon1,recMuon2);
+           AliAODDimuon *dimuonJpsi = new AliAODDimuon;
+           dimuonJpsi->SetMuons(recMuon1,recMuon2);
            Double_t pca[3]={0};
            Double_t pcaQuality=0;
            TLorentzVector kinem(0,0,0,0);
-           if (!AliMFTAnalysisTools::CalculatePCA(dimuon, pca, pcaQuality, kinem)) continue;
+           if (!AliMFTAnalysisTools::CalculatePCA(dimuonJpsi, pca, pcaQuality, kinem)) {
+             delete dimuonJpsi;
+             continue;
+           }
            fHistPtDimuonsJpsi    -> Fill(kinem.Pt());
            fHistMassDimuonsJpsi  -> Fill(kinem.M());
            fHistResidualXVtxJpsi -> Fill(1.e4*(pca[0] - fPrimaryVertex[0]));
            fHistResidualYVtxJpsi -> Fill(1.e4*(pca[1] - fPrimaryVertex[1]));
            fHistResidualZVtxJpsi -> Fill(1.e4*(pca[2] - fPrimaryVertex[2]));
-           delete dimuon;
+           delete dimuonJpsi;
          }
        }
       }
index 815d929646b824352038d6e662c6e07c4b17eed3..d65930879fbd6c4891ae717ab87adfde716e9cb2 100644 (file)
@@ -137,8 +137,6 @@ Bool_t AliMFTAnalysisTools::ExtrapAODMuonToXY(AliAODTrack *muon, Double_t xy[2],
 
   // We look for the above-defined PCA
 
-  Double_t xPCA=0, yPCA=0, zPCA=0;
-
   AliMUONTrackParam *param = new AliMUONTrackParam();
   param -> SetNonBendingCoor(muon->XAtDCA());
   param -> SetBendingCoor(muon->YAtDCA());
@@ -675,3 +673,80 @@ const TMatrixD AliMFTAnalysisTools::ConvertCovMatrixAOD2MUON(AliAODTrack *muon)
 
 //====================================================================================================================================================
 
+Bool_t AliMFTAnalysisTools::IsCorrectMatch(AliAODTrack *muon) {
+
+  for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) if (IsWrongCluster(muon, iPlane)) return kFALSE;
+  return kTRUE;
+
+}
+
+//====================================================================================================================================================
+
+TString AliMFTAnalysisTools::GetGenerator(Int_t label, AliAODMCHeader* header) {
+
+  // get the name of the generator that produced a given particle
+
+  Int_t partCounter = 0;
+  TList *genHeaders = header->GetCocktailHeaders();
+  Int_t nGenHeaders = genHeaders->GetEntries();
+
+  for (Int_t i=0; i<nGenHeaders; i++){
+    AliGenEventHeader *gh = (AliGenEventHeader*) genHeaders->At(i);
+    TString genName = gh->GetName();
+    Int_t nPart = gh->NProduced();
+    if (label>=partCounter && label<(partCounter+nPart)) return genName;
+    partCounter += nPart;
+  }
+
+  TString empty="";
+  return empty;
+
+}
+
+//====================================================================================================================================================
+
+void AliMFTAnalysisTools::GetTrackPrimaryGenerator(AliAODTrack *track, AliAODMCHeader *header, TClonesArray *arrayMC, TString &nameGen) {
+
+  // method to check if a track comes from a given generator
+
+  Int_t label = TMath::Abs(track->GetLabel());
+  nameGen = GetGenerator(label,header);
+  
+  // In case the particle is not primary nameGen will contain blank spaces. In this case, we search backward for the primary which originated the chain
+  
+  while (nameGen.IsWhitespace()) {
+    AliAODMCParticle *mcPart = (AliAODMCParticle*) arrayMC->At(label);
+    if (!mcPart) {
+      printf("AliMFTAnalysisTools::GetTrackPrimaryGenerator - BREAK: No valid AliAODMCParticle at label %i\n",label);
+      break;
+    }
+    Int_t motherLabel = mcPart->GetMother();
+    if (motherLabel < 0) {
+      printf("AliMFTAnalysisTools::GetTrackPrimaryGenerator - BREAK: Reached primary particle without valid mother\n");
+      break;
+    }
+    label = motherLabel;
+    nameGen = GetGenerator(label,header);
+  }
+  
+  return;
+
+}
+
+//====================================================================================================================================================
+
+Bool_t AliMFTAnalysisTools::IsTrackInjected(AliAODTrack *track, AliAODMCHeader *header, TClonesArray *arrayMC){
+
+  // method to check if a track comes from the signal event or from the underlying Hijing event
+
+  TString nameGen;
+
+  GetTrackPrimaryGenerator(track,header,arrayMC,nameGen);
+  
+  if (nameGen.IsWhitespace() || nameGen.Contains("ijing")) return kFALSE;
+  
+  return kTRUE;
+}
+
+//====================================================================================================================================================
+
index b51a633956e3566cf97f920567f9d9bc521915db..e6e0a7f7c2ef0f7737088b905152f1d6ebc2a47f 100644 (file)
@@ -22,6 +22,9 @@
 #include "AliLog.h"
 #include "TMatrixD.h"
 #include "TClonesArray.h"
+#include "AliAODMCHeader.h"
+#include "AliGenEventHeader.h"
+#include "AliAODMCParticle.h"
 
 //====================================================================================================================================================
 
@@ -58,6 +61,12 @@ public:
     else return !(muon->GetMFTClusterPattern() & (1<<(iPlane+AliMFTConstants::fNMaxPlanes)));
   }
 
+  static Bool_t IsCorrectMatch(AliAODTrack *muon);
+
+  static TString GetGenerator(Int_t label, AliAODMCHeader* header);
+  static void GetTrackPrimaryGenerator(AliAODTrack *track, AliAODMCHeader *header, TClonesArray *arrayMC, TString &nameGen);
+  static Bool_t IsTrackInjected(AliAODTrack *track, AliAODMCHeader *header, TClonesArray *arrayMC);
+
   static void ConvertCovMatrixMUON2AOD(const TMatrixD& covMUON, Double_t covAOD[21]);
   static const TMatrixD ConvertCovMatrixAOD2MUON(AliAODTrack *muon);
   
index bad7caa944dbab6b194bbd0cdd9e5d3afd250bd9..9a5dbd94f93548b5920428a8340bb08832f8032a 100644 (file)
@@ -683,7 +683,7 @@ void AliMFTPlane::DrawPlane(Option_t *opt) {
     cnv->Draw();
 
     TH2D *h = new TH2D("tmp", GetName(), 
-                      1, fZCenter-0.5, fZCenter+0.5, 
+                      1, fZCenter-1.1*(0.5*fThicknessSupport+fThicknessActive), fZCenter+1.1*(0.5*fThicknessSupport+fThicknessActive),
                       1, 1.1*GetSupportElement(0)->GetAxis(1)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(1)->GetXmax());
     h->SetXTitle("z [cm]");
     h->SetYTitle("y [cm]");