]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Possibility to vary the cut on DCA between V0 daughters
authorbelikov <belikov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 23 Mar 2012 17:30:32 +0000 (17:30 +0000)
committerbelikov <belikov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 23 Mar 2012 17:30:32 +0000 (17:30 +0000)
PWGLF/STRANGENESS/LambdaK0PbPb/AddTaskCTauAOD.C
PWGLF/STRANGENESS/LambdaK0PbPb/AliAnalysisTaskCTauPbPbaod.cxx
PWGLF/STRANGENESS/LambdaK0PbPb/AliAnalysisTaskCTauPbPbaod.h

index b9e0cd53d09fefee5e8d1815d77196833a34dc15..ac0b7f38092c0b7c8f8afb41a897d98a0dda7bf3 100644 (file)
@@ -1,6 +1,6 @@
 AliAnalysisTaskCTauPbPbaod* 
 AddTaskCTauAOD(Double_t min=0., Double_t max=90., Double_t cpa=0.9975, 
-TString name="cTau_0090aod", Bool_t isMC=kFALSE) 
+Double_t dca=1.5, TString name="cTau_0090aod", Bool_t isMC=kFALSE) 
 {
   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
   if (!mgr) {
@@ -14,8 +14,9 @@ TString name="cTau_0090aod", Bool_t isMC=kFALSE)
   }
 
   AliAnalysisTaskCTauPbPbaod *task = new AliAnalysisTaskCTauPbPbaod(name);
-  task->SetCosPA(cpa);
   task->SetCentrality(min,max);
+  task->SetCosPA(cpa);
+  task->SetDtrDCA(dca);
   task->SetMC(isMC);
   mgr->AddTask(task);
   
index a1d75239f8db6c6974a69fb50fcc5846ae450b29..680bcb31a50405096d4584b01f3d4bda16818640 100644 (file)
@@ -42,9 +42,11 @@ fIsMC(kFALSE),
 fCMin(0.),
 fCMax(90.),
 fCPA(0.9975),
+fDCA(1.0),
 fOutput(0),
 fMult(0),
 fCosPA(0),
+fDtrDCA(0),
 fdEdx(0),
 fdEdxPid(0),
 
@@ -99,6 +101,10 @@ void AliAnalysisTaskCTauPbPbaod::UserCreateOutputObjects()
   fCosPA->GetXaxis()->SetTitle("Cos(PA)"); 
   fOutput->Add(fCosPA);
 
+  fDtrDCA=new TH1F("fDtrDCA","DCA between V0 daughters",50,0.0,1.5);
+  fDtrDCA->GetXaxis()->SetTitle("DCA (rel. u.)"); 
+  fOutput->Add(fDtrDCA);
+
   fdEdx=new TH2F("fdEdx","dE/dx",50,0.2,3,50,0.,6.);
   fOutput->Add(fdEdx);
 
@@ -266,6 +272,9 @@ AliAnalysisTaskCTauPbPbaod::AcceptV0(const AliAODv0 *v0,const AliAODEvent *aod)
   Double_t cpa=v0->CosPointingAngle(aod->GetPrimaryVertex());
   if (cpa < fCPA) return kFALSE;
 
+  Double_t dca=v0->DcaV0Daughters();
+  if (dca > fDCA) return kFALSE;
+
   const AliAODTrack *ntrack=(AliAODTrack *)v0->GetDaughter(1);
   if (!AcceptTrack(ntrack)) return kFALSE;
 
@@ -277,11 +286,6 @@ AliAnalysisTaskCTauPbPbaod::AcceptV0(const AliAODv0 *v0,const AliAODEvent *aod)
   xy=v0->DcaPosToPrimVertex();
   if (TMath::Abs(xy)<0.1) return kFALSE;
 
-  Double_t dca=v0->DcaV0Daughters();
-  if (dca>1.0) return kFALSE;
-  //if (dca>0.7) return kFALSE;
-  //if (dca>0.4) return kFALSE;
-
   Double_t xyz[3]; v0->GetSecondaryVtx(xyz);
   Double_t r2=xyz[0]*xyz[0] + xyz[1]*xyz[1];
   if (r2<0.9*0.9) return kFALSE;
@@ -515,6 +519,9 @@ void AliAnalysisTaskCTauPbPbaod::UserExec(Option_t *)
       Double_t cpa=v0->CosPointingAngle(aod->GetPrimaryVertex());
       fCosPA->Fill(cpa);
 
+      Double_t dca=v0->DcaV0Daughters();
+      fDtrDCA->Fill(dca);
+
       const AliAODTrack *ntrack=(AliAODTrack *)v0->GetDaughter(1);
       const AliAODTrack *ptrack=(AliAODTrack *)v0->GetDaughter(0);
 
index 9ee02fc6e2713da15a94e85bccbaba9f4914be72..f4fb90971c37c9dd8d60e474fcf9c71a2d493428 100644 (file)
@@ -25,6 +25,7 @@ public:
   void SetCentrality(Double_t min, Double_t max) {fCMin=min;fCMax=max;} 
   void SetMC(Bool_t isMC=kTRUE) {fIsMC=isMC;} 
   void SetCosPA(Double_t cospa) {fCPA=cospa;} 
+  void SetDtrDCA(Double_t cospa){fDCA=cospa;} 
   
   virtual void   UserCreateOutputObjects();
   virtual void   UserExec(Option_t *option);
@@ -42,11 +43,13 @@ private:
   Double_t fCMin;       // Min centrality
   Double_t fCMax;       // Max centrality
   Double_t fCPA;        // cos(PA) threshold
+  Double_t fDCA;        // threshold for the DCA between V0 daughters
 
   TList       *fOutput; //! The list of histograms
 
   TH1F *fMult;       //! Track multiplicity
-  TH1F *fCosPA;      //! Track multiplicity
+  TH1F *fCosPA;      //! cos(PA)
+  TH1F *fDtrDCA;     //! DCA between V0 daughters
   TH2F* fdEdx;       //! dEdx
   TH2F* fdEdxPid;    //! dEdx with PID
 
@@ -81,7 +84,7 @@ private:
   TH2F* fXiBarM;         //! Mass for anti-Xis
   TH1F* fXiBarSiP;       //! Side-band subtracted Pt for reconstructed anti-Xi
 
-  ClassDef(AliAnalysisTaskCTauPbPbaod,4);
+  ClassDef(AliAnalysisTaskCTauPbPbaod,5);
 };
 
 #endif