]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
update aod track cuts
authorSvein Lindal <svein.lindal@gmail.com>
Sun, 8 Dec 2013 23:00:05 +0000 (18:00 -0500)
committerSvein Lindal <svein.lindal@gmail.com>
Sun, 8 Dec 2013 23:07:01 +0000 (18:07 -0500)
change track id checking for aods
modified:   AliAnalysisTaskdPhi.cxx
modified:   AliConversionTrackCuts.cxx
modified:   AliConversionTrackCuts.h

PWGGA/GammaConv/AliAnalysisTaskdPhi.cxx
PWGGA/GammaConv/AliConversionTrackCuts.cxx
PWGGA/GammaConv/AliConversionTrackCuts.h

index 2642503edcdf395b3202519593285a5a453e198e..1414546e34a43c44cde9a8bfc59325e7ac4a9563 100644 (file)
@@ -689,8 +689,6 @@ void AliAnalysisTaskdPhi::UserExec(Option_t *) {
   Bool_t lpitmap[fAxistPt.GetNbins()][fMesonFilters[0].GetEntriesFast()];
   Bool_t upitmap[fAxistPt.GetNbins()][fMesonFilters[1].GetEntriesFast()];
   
-  
-
   for(Int_t igf = 0; igf < fV0Filters[0].GetEntriesFast(); igf++) {
     for(Int_t ptbin = 0; ptbin < fAxistPt.GetNbins(); ptbin++) {
       lv0tmap[ptbin][igf] = kFALSE;
@@ -899,10 +897,18 @@ void AliAnalysisTaskdPhi::UserExec(Option_t *) {
                  for(int ij = 0; ij < tracks->GetEntriesFast(); ij++) {
                    AliVTrack * track = static_cast<AliVTrack*>(tracks->At(ij));
                    Int_t tid = track->GetID();
-                 
-                   if((tid > 0) && (tid == tIDs[0] || tid == tIDs[1] || tid == tIDs[2] || tid == tIDs[3]) ) {
+
+
+                   if(tid == tIDs[0] || tid == tIDs[1] || tid == tIDs[2] || tid == tIDs[3])  {
                      continue;
                    }
+                   
+                   if(tid < 0) {
+                     if(-tid-1 == tIDs[0]+1 || -tid-1 == tIDs[1]+1 || -tid-1 == tIDs[2]+1 || -tid-1 == tIDs[3]+1)  {
+                       continue;
+                     }
+                   }
+                 
                  
                    dphivalues[0] = pion->Eta() - track->Eta();
                    dphivalues[1] = GetDPhi(pion->Phi() - track->Phi());
@@ -980,10 +986,16 @@ void AliAnalysisTaskdPhi::UserExec(Option_t *) {
                        AliVTrack * track = static_cast<AliVTrack*>(ttracks->At(ij));
                        Int_t tid = track->GetID();
                        
-                       if((tid > 0) && (tid == tIDs[0] || tid == tIDs[1] || tid == tIDs[2] || tid == tIDs[3]) ) {
+                       if(tid == tIDs[0] || tid == tIDs[1] || tid == tIDs[2] || tid == tIDs[3] ) {
                          continue;
                        }
-                       
+
+                       if(tid < 0) {
+                         if(-tid-1 == tIDs[0]+1 || -tid-1 == tIDs[1]+1 || -tid-1 == tIDs[2]+1 || -tid-1 == tIDs[3]+1)  {
+                           continue;
+                         }
+                       }
+                 
                        dphivalues[0] = pion->Eta() - track->Eta();
                        dphivalues[1] = GetDPhi(pion->Phi() - track->Phi());
                        dphivalues[3] = track->Pt();
@@ -1068,7 +1080,7 @@ void AliAnalysisTaskdPhi::UserExec(Option_t *) {
       for(Int_t iguf = 0; iguf < fV0Filters[1].GetEntriesFast(); iguf++) {
        if (uv0tmap[tbin][iguf] ) {
 
-         cout << "c vtx " <<  centrality << vertexz << endl;
+         //cout << "c vtx " <<  centrality << vertexz << endl;
 
          for(Int_t iTrack = 0; iTrack < ttracks->GetEntriesFast(); iTrack++) {
            AliVTrack * track = static_cast<AliVTrack*>(ttracks->At(iTrack));
@@ -1085,7 +1097,7 @@ void AliAnalysisTaskdPhi::UserExec(Option_t *) {
       for(Int_t ipuf = 0; ipuf < fMesonFilters[1].GetEntriesFast(); ipuf++) {
        if (upitmap[tbin][ipuf] ) {
 
-         cout << "c2 vtx2 " <<  centrality << vertexz << endl;
+         //cout << "c2 vtx2 " <<  centrality << vertexz << endl;
 
 
          for(Int_t iTrack = 0; iTrack < ttracks->GetEntriesFast(); iTrack++) {
@@ -1099,7 +1111,6 @@ void AliAnalysisTaskdPhi::UserExec(Option_t *) {
          }
        }
       }
-
     }
   }
 
index fa0495009f866b946f76e3d7308a3b3a92e857bc..47a6ef5136535ead03d880d344411cb67bca5bd6 100644 (file)
@@ -56,13 +56,14 @@ AliConversionTrackCuts::AliConversionTrackCuts() :
   fOwnedTracks(),
   fInitialized(kFALSE),
   fhPhi(NULL),
-  fhPt(NULL),
-  fhPhiPt(NULL),
+  //  fhPt(NULL),
+  //fhPhiPt(NULL),
   fhdcaxyPt(NULL),
   fhdcazPt(NULL),
   fhdca(NULL),
   fhnclpt(NULL),
   fhnclsfpt(NULL),
+  fhEtaPhi(NULL),
   fHistograms(NULL) 
 {
   //Constructor
@@ -81,13 +82,14 @@ AliConversionTrackCuts::AliConversionTrackCuts(TString name, TString title = "ti
   fOwnedTracks(),
   fInitialized(kFALSE),
   fhPhi(NULL),  
-  fhPt(NULL),
-  fhPhiPt(NULL),
+  //fhPt(NULL),
+  //fhPhiPt(NULL),
   fhdcaxyPt(NULL),
   fhdcazPt(NULL),
   fhdca(NULL),
   fhnclpt(NULL),
   fhnclsfpt(NULL),
+  fhEtaPhi(NULL),
   fHistograms(NULL)
 {
   //Constructor
@@ -275,57 +277,77 @@ Bool_t AliConversionTrackCuts::AcceptTrack(AliESDtrack * track) {
 
 Bool_t AliConversionTrackCuts::AcceptTrack(AliAODTrack * track) {
   //Check aod track
+
   FillHistograms(kPreCut, track);
-  if(!track->TestFilterBit(fFilterBit)) {
-    return kFALSE;
-  }
 
+  if(!track->IsHybridGlobalConstrainedGlobal()) return kFALSE;
 
+   if (!(track->GetStatus() & AliVTrack::kITSrefit)) {
+     return kFALSE;
+   }
 
 
-  ///Do dca xy cut!
-  FillHistograms(1, track);
-  return kTRUE;
+  //The cluster sharing cut can be done with:
+  Double_t frac = Double_t(track->GetTPCnclsS()) / Double_t(track->GetTPCncls());
+  if (frac > 0.4) return kFALSE;
 
 
-  // if (track->GetTPCNcls() < fTPCminNClusters) return kFALSE;
-  // FillHistograms(kCutNcls, track);
-
-  // if (track->Chi2perNDF() > fTPCmaxChi2) return kFALSE;
-  // FillHistograms(kCutNDF, track);
+  ///Do dca xy cut!
+  FillHistograms(1, track);
 
-  // AliAODVertex *vertex = track->GetProdVertex();
-  // if (vertex && fRejectKinkDaughters) {
-  //   if (vertex->GetType() == AliAODVertex::kKink) {
-  //     return kFALSE;
-  //   }
-  // }
-  // FillHistograms(kCutKinc, track);
+  ///DCA
+  Double_t dca[2] = { -999, -999};
+  //Bool_t dcaok = 
+  GetDCA(track, dca);
+  FillDCAHist(dca[1], dca[0], track);
   
-  // if(TMath::Abs(track->ZAtDCA()) > fDCAZmax) {
-  //   return kFALSE;
-  // }
-  // FillHistograms(kCutDCAZ, track);
-
-  // Float_t xatdca = track->XAtDCA();
-  // Float_t yatdca = track->YAtDCA();
-  // Float_t xy = xatdca*xatdca + yatdca*yatdca;
-  // if(xy > fDCAXYmax) {
-  //   return kFALSE;
-  // }
-
-  // FillHistograms(kCutDCAXY, track);
-
-
-  // fhnclpt->Fill(track->Pt(), track->GetTPCNcls());
-  // if(track->GetTPCNclsF() > 0) fhnclsfpt->Fill(track->Pt(), ((Double_t) track->GetTPCNcls())/track->GetTPCNclsF());
-  // FillDCAHist(track->ZAtDCA(), TMath::Sqrt(track->XAtDCA()*track->XAtDCA() + track->YAtDCA()*track->YAtDCA()), track);
 
+  if(track->IsGlobalConstrained()) {
+    FillHistograms(3, track);
+  } else {
+    FillHistograms(2, track);
+  }
+  
+  if(fhEtaPhi) fhEtaPhi->Fill(track->Eta(), track->Phi());
 
+  return kTRUE;
 }
 
+///______________________________________________________________________________
+Bool_t AliConversionTrackCuts::GetDCA(const AliAODTrack *track, Double_t dca[2]) {
+  if(track->TestBit(AliAODTrack::kIsDCA)){
+    dca[0]=track->DCA();
+    dca[1]=track->ZAtDCA();
+    return kTRUE;
+  }
+  
+  Bool_t ok=kFALSE;
+  if(fEvent) {
+    Double_t covdca[3];
+    //AliAODTrack copy(*track);
+    AliExternalTrackParam etp; etp.CopyFromVTrack(track);
+    
+    Float_t xstart = etp.GetX();
+    if(xstart>3.) {
+      dca[0]=-999.;
+    dca[1]=-999.;
+    //printf("This method can be used only for propagation inside the beam pipe \n");
+    return kFALSE;
+    }
 
 
+    AliAODVertex *vtx =(AliAODVertex*)(fEvent->GetPrimaryVertex());
+    Double_t fBzkG = fEvent->GetMagneticField(); // z componenent of field in kG
+    ok = etp.PropagateToDCA(vtx,fBzkG,kVeryBig,dca,covdca);
+    //ok = copy.PropagateToDCA(vtx,fBzkG,kVeryBig,dca,covdca);
+  }
+  if(!ok){
+    dca[0]=-999.;
+    dca[1]=-999.;
+  }
+  return ok;
+}
 
 
 
@@ -345,6 +367,8 @@ TList * AliConversionTrackCuts::CreateHistograms() {
   fHistograms->Add(fhPhi);
   
 
+  fhEtaPhi = new TH2F("etahpi", "etaphi", 36, -0.9, 0.9, 32, 0, TMath::TwoPi());
+  fHistograms->Add(fhEtaPhi);
 
   // fhPt = new TH2F("pt", "pt", kNCuts+2, kPreCut -0.5, kNCuts + 0.5, 
   //                             20, 0., 20.);
@@ -357,14 +381,18 @@ TList * AliConversionTrackCuts::CreateHistograms() {
   //  fhPhiPt = new TH2F("phipt", "phipt", 100, 0, 100, 64, 0, TMath::TwoPi());
   //fHistograms->Add(fhPhiPt);
 
-  // fhdcaxyPt = new TH2F("dcaxypt", "dcaxypt", 20, 0, 20, 50, 0, 5);
-  // fHistograms->Add(fhdcaxyPt);
+  fhdcaxyPt = new TH2F("dcaxypt", "dcaxypt", 20, 0, 20, 50, -2.5, 2.5);
+  fHistograms->Add(fhdcaxyPt);
 
-  // fhdcazPt = new TH2F("dcazpt", "dcazpt", 20, 0, 20, 50, 0, 5);
-  // fHistograms->Add(fhdcazPt);
+  fhdcazPt = new TH2F("dcazpt", "dcazpt", 20, 0, 20, 70, -3.5, 3.5);
+  fHistograms->Add(fhdcazPt);
 
-  // fhdca = new TH2F("dca", "dca", 60, -3, 3, 60, -3, 3);
-  // fHistograms->Add(fhdca);
+  fhdca = new TH2F("dca", "dca", 70, -3.5, 3.5, 50, -2.5, 2.5);
+  fhdca->SetXTitle("dca z");
+  fhdca->SetYTitle("dca xy");
+
+  
+  fHistograms->Add(fhdca);
 
   // fhnclpt = new TH2F("nclstpcvspt", "nclstpcvspt", 20, 0, 20, 50, 0, 100);
   // fHistograms->Add(fhnclpt);
index 195cd00781c4954c29bd86d67f8df555b202a774..03a3f930eceae7f677b76f8e3f4e805d93e98b7c 100644 (file)
@@ -42,6 +42,8 @@ public:
   Bool_t IsSelected(TList * /*list*/) { return kFALSE; }
   Bool_t AcceptTrack(AliAODTrack * track);
   Bool_t AcceptTrack(AliESDtrack * track);
+  Bool_t GetDCA(const AliAODTrack * track, Double_t dca[2]);
+
 
   void DeleteTracks() { fOwnedTracks.Delete(); }
   void FillDCAHist(Float_t dcaz, Float_t dcaxy, AliVTrack * track);
@@ -82,20 +84,21 @@ protected :
   Bool_t fInitialized;
 
   TH2F * fhPhi;
-  TH2F * fhPt;
-  TH2F * fhPhiPt;
+  //TH2F * fhPt;
+  //TH2F * fhPhiPt;
   TH2F * fhdcaxyPt;
   TH2F * fhdcazPt;
   TH2F * fhdca;
   TH2F * fhnclpt;
   TH2F * fhnclsfpt;
+  TH2F * fhEtaPhi;
   
   TList * fHistograms;
 
   AliConversionTrackCuts(const AliConversionTrackCuts&); // not implemented
   AliConversionTrackCuts& operator=(const AliConversionTrackCuts&); // not implemented
 
-  ClassDef(AliConversionTrackCuts,3)
+  ClassDef(AliConversionTrackCuts,4)
 };
 
 #endif