]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Lightweight cascade checks -> code implementation 2
authorddobrigk <david.dobrigkeit.chinellato@cern.ch>
Wed, 26 Nov 2014 19:09:48 +0000 (17:09 -0200)
committerddobrigk <david.dobrigkeit.chinellato@cern.ch>
Wed, 26 Nov 2014 19:10:36 +0000 (17:10 -0200)
PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightCascadeVertexer.cxx [new file with mode: 0644]
PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightCascadeVertexer.h [new file with mode: 0644]
PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightV0vertexer.cxx [new file with mode: 0644]
PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightV0vertexer.h [new file with mode: 0644]

diff --git a/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightCascadeVertexer.cxx b/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightCascadeVertexer.cxx
new file mode 100644 (file)
index 0000000..a4bd5f3
--- /dev/null
@@ -0,0 +1,271 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+//-------------------------------------------------------------------------
+//               Implementation of the cascade vertexer class
+//          Reads V0s and tracks, writes out cascade vertices
+//                     Fills the ESD with the cascades 
+//    Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
+//-------------------------------------------------------------------------
+
+//modified by R. Vernet 30/6/2006 : daughter label
+//modified by R. Vernet  3/7/2006 : causality
+//modified by I. Belikov 24/11/2006 : static setter for the default cuts
+
+#include "AliESDEvent.h"
+#include "AliESDcascade.h"
+#include "AliLightCascadeVertexer.h"
+
+ClassImp(AliLightCascadeVertexer)
+
+//A set of loose cuts
+Double_t 
+  AliLightCascadeVertexer::fgChi2max=33.;   //maximal allowed chi2 
+Double_t 
+  AliLightCascadeVertexer::fgDV0min=0.01;   //min V0 impact parameter
+Double_t 
+  AliLightCascadeVertexer::fgMassWin=0.008; //"window" around the Lambda mass
+Double_t 
+  AliLightCascadeVertexer::fgDBachMin=0.01; //min bachelor impact parameter
+Double_t 
+  AliLightCascadeVertexer::fgDCAmax=2.0;    //max DCA between the V0 and the track 
+Double_t 
+  AliLightCascadeVertexer::fgCPAmin=0.98; //min cosine of the cascade pointing angle
+Double_t 
+  AliLightCascadeVertexer::fgRmin=0.2;      //min radius of the fiducial volume
+Double_t 
+  AliLightCascadeVertexer::fgRmax=100.;     //max radius of the fiducial volume
+
+Double_t AliLightCascadeVertexer::fgMaxEta=0.8;        //max |eta|
+Double_t AliLightCascadeVertexer::fgMinClusters=70;   //min clusters (>=)
+  
+
+Int_t AliLightCascadeVertexer::V0sTracks2CascadeVertices(AliESDEvent *event) {
+  //--------------------------------------------------------------------
+  // This function reconstructs cascade vertices
+  //      Adapted to the ESD by I.Belikov (Jouri.Belikov@cern.ch)
+  //--------------------------------------------------------------------
+   const AliESDVertex *vtxT3D=event->GetPrimaryVertex();
+
+   Double_t xPrimaryVertex=vtxT3D->GetX();
+   Double_t yPrimaryVertex=vtxT3D->GetY();
+   Double_t zPrimaryVertex=vtxT3D->GetZ();
+
+   Double_t b=event->GetMagneticField();
+   Int_t nV0=(Int_t)event->GetNumberOfV0s();
+
+   //stores relevant V0s in an array
+   TObjArray vtcs(nV0);
+   Int_t i;
+   for (i=0; i<nV0; i++) {
+       AliESDv0 *v=event->GetV0(i);
+       if (v->GetOnFlyStatus()) continue;
+       if (v->GetD(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex)<fDV0min) continue;
+       vtcs.AddLast(v);
+   }
+   nV0=vtcs.GetEntriesFast();
+
+   // stores relevant tracks in another array
+   Int_t nentr=(Int_t)event->GetNumberOfTracks();
+   TArrayI trk(nentr); Int_t ntr=0;
+   for (i=0; i<nentr; i++) {
+       AliESDtrack *esdtr=event->GetTrack(i);
+       ULong_t status=esdtr->GetStatus();
+
+       if ((status&AliESDtrack::kITSrefit)==0)
+         if ((status&AliESDtrack::kTPCrefit)==0) continue;
+       
+       //Track pre-selection: clusters
+       if (esdtr->GetTPCNcls() < fMinClusters ) continue;
+
+       if (TMath::Abs(esdtr->GetD(xPrimaryVertex,yPrimaryVertex,b))<fDBachMin) continue;
+
+       trk[ntr++]=i;
+   }   
+
+   Double_t massLambda=1.11568;
+   Int_t ncasc=0;
+
+   // Looking for the cascades...
+
+   for (i=0; i<nV0; i++) { //loop on V0s
+
+      AliESDv0 *v=(AliESDv0*)vtcs.UncheckedAt(i);
+      AliESDv0 v0(*v);
+      v0.ChangeMassHypothesis(kLambda0); // the v0 must be Lambda 
+      if (TMath::Abs(v0.GetEffMass()-massLambda)>fMassWin) continue; 
+
+      for (Int_t j=0; j<ntr; j++) {//loop on tracks
+        Int_t bidx=trk[j];
+        //Bo:   if (bidx==v->GetNindex()) continue; //bachelor and v0's negative tracks must be different
+         if (bidx==v0.GetIndex(0)) continue; //Bo:  consistency 0 for neg
+        AliESDtrack *btrk=event->GetTrack(bidx);
+         if (btrk->GetSign()>0) continue;  // bachelor's charge 
+          
+        AliESDv0 *pv0=&v0;
+         AliExternalTrackParam bt(*btrk), *pbt=&bt;
+
+         Double_t dca=PropagateToDCA(pv0,pbt,b);
+         if (dca > fDCAmax) continue;
+          
+          //eta cut - test
+            if (TMath::Abs(pbt->Eta())>fMaxEta) continue;
+
+         AliESDcascade cascade(*pv0,*pbt,bidx);//constucts a cascade candidate
+        //PH        if (cascade.GetChi2Xi() > fChi2max) continue;
+
+        Double_t x,y,z; cascade.GetXYZcascade(x,y,z); // Bo: bug correction
+         Double_t r2=x*x + y*y; 
+         if (r2 > fRmax*fRmax) continue;   // condition on fiducial zone
+         if (r2 < fRmin*fRmin) continue;
+
+        Double_t pxV0,pyV0,pzV0;
+        pv0->GetPxPyPz(pxV0,pyV0,pzV0);
+        if (x*pxV0+y*pyV0+z*pzV0 < 0) continue; //causality
+
+         Double_t x1,y1,z1; pv0->GetXYZ(x1,y1,z1);
+         if (r2 > (x1*x1+y1*y1)) continue;
+
+        if (cascade.GetCascadeCosineOfPointingAngle(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex) <fCPAmin) continue; //condition on the cascade pointing angle 
+        
+         cascade.SetDcaXiDaughters(dca);
+        event->AddCascade(&cascade);
+         ncasc++;
+      } // end loop tracks
+   } // end loop V0s
+
+   // Looking for the anti-cascades...
+
+   for (i=0; i<nV0; i++) { //loop on V0s
+      AliESDv0 *v=(AliESDv0*)vtcs.UncheckedAt(i);
+      AliESDv0 v0(*v);
+      v0.ChangeMassHypothesis(kLambda0Bar); //the v0 must be anti-Lambda 
+      if (TMath::Abs(v0.GetEffMass()-massLambda)>fMassWin) continue; 
+
+      for (Int_t j=0; j<ntr; j++) {//loop on tracks
+        Int_t bidx=trk[j];
+        //Bo:   if (bidx==v->GetPindex()) continue; //bachelor and v0's positive tracks must be different
+         if (bidx==v0.GetIndex(1)) continue; //Bo:  consistency 1 for pos
+        AliESDtrack *btrk=event->GetTrack(bidx);
+         if (btrk->GetSign()<0) continue;  // bachelor's charge 
+          
+        AliESDv0 *pv0=&v0;
+         AliExternalTrackParam bt(*btrk), *pbt=&bt;
+
+         Double_t dca=PropagateToDCA(pv0,pbt,b);
+         if (dca > fDCAmax) continue;
+
+         AliESDcascade cascade(*pv0,*pbt,bidx); //constucts a cascade candidate
+        //PH         if (cascade.GetChi2Xi() > fChi2max) continue;
+
+        Double_t x,y,z; cascade.GetXYZcascade(x,y,z); // Bo: bug correction
+         Double_t r2=x*x + y*y; 
+         if (r2 > fRmax*fRmax) continue;   // condition on fiducial zone
+         if (r2 < fRmin*fRmin) continue;
+
+        Double_t pxV0,pyV0,pzV0;
+        pv0->GetPxPyPz(pxV0,pyV0,pzV0);
+        if (x*pxV0+y*pyV0+z*pzV0 < 0) continue; //causality
+
+         Double_t x1,y1,z1; pv0->GetXYZ(x1,y1,z1);
+         if (r2 > (x1*x1+y1*y1)) continue;
+
+        if (cascade.GetCascadeCosineOfPointingAngle(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex) < fCPAmin) continue; //condition on the cascade pointing angle 
+
+         cascade.SetDcaXiDaughters(dca);
+        event->AddCascade(&cascade);
+         ncasc++;
+
+      } // end loop tracks
+   } // end loop V0s
+
+Info("V0sTracks2CascadeVertices","Number of reconstructed cascades: %d",ncasc);
+
+   return 0;
+}
+
+
+Double_t AliLightCascadeVertexer::Det(Double_t a00, Double_t a01, Double_t a10, Double_t a11) const {
+  //--------------------------------------------------------------------
+  // This function calculates locally a 2x2 determinant
+  //--------------------------------------------------------------------
+  return a00*a11 - a01*a10;
+}
+
+Double_t AliLightCascadeVertexer::Det(Double_t a00,Double_t a01,Double_t a02,
+                                Double_t a10,Double_t a11,Double_t a12,
+                                Double_t a20,Double_t a21,Double_t a22) const {
+  //--------------------------------------------------------------------
+  // This function calculates locally a 3x3 determinant
+  //--------------------------------------------------------------------
+  return  a00*Det(a11,a12,a21,a22)-a01*Det(a10,a12,a20,a22)+a02*Det(a10,a11,a20,a21);
+}
+
+
+
+
+Double_t AliLightCascadeVertexer::PropagateToDCA(AliESDv0 *v, AliExternalTrackParam *t, Double_t b) {
+  //--------------------------------------------------------------------
+  // This function returns the DCA between the V0 and the track
+  //--------------------------------------------------------------------
+  Double_t alpha=t->GetAlpha(), cs1=TMath::Cos(alpha), sn1=TMath::Sin(alpha);
+  Double_t r[3]; t->GetXYZ(r);
+  Double_t x1=r[0], y1=r[1], z1=r[2];
+  Double_t p[3]; t->GetPxPyPz(p);
+  Double_t px1=p[0], py1=p[1], pz1=p[2];
+  
+  Double_t x2,y2,z2;     // position and momentum of V0
+  Double_t px2,py2,pz2;
+  
+  v->GetXYZ(x2,y2,z2);
+  v->GetPxPyPz(px2,py2,pz2);
+// calculation dca
+   
+  Double_t dd= Det(x2-x1,y2-y1,z2-z1,px1,py1,pz1,px2,py2,pz2);
+  Double_t ax= Det(py1,pz1,py2,pz2);
+  Double_t ay=-Det(px1,pz1,px2,pz2);
+  Double_t az= Det(px1,py1,px2,py2);
+
+  Double_t dca=TMath::Abs(dd)/TMath::Sqrt(ax*ax + ay*ay + az*az);
+
+//points of the DCA
+  Double_t t1 = Det(x2-x1,y2-y1,z2-z1,px2,py2,pz2,ax,ay,az)/
+                Det(px1,py1,pz1,px2,py2,pz2,ax,ay,az);
+  
+  x1 += px1*t1; y1 += py1*t1; //z1 += pz1*t1;
+  
+
+  //propagate track to the points of DCA
+
+  x1=x1*cs1 + y1*sn1;
+  if (!t->PropagateTo(x1,b)) {
+    Error("PropagateToDCA","Propagation failed !");
+    return 1.e+33;
+  }  
+
+  return dca;
+}
+
+
+
+
+
+
+
+
+
+
+
diff --git a/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightCascadeVertexer.h b/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightCascadeVertexer.h
new file mode 100644 (file)
index 0000000..50fd358
--- /dev/null
@@ -0,0 +1,129 @@
+#ifndef AliLightCascadeVertexer_H
+#define AliLightCascadeVertexer_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//------------------------------------------------------------------
+//                    Cascade Vertexer Class
+//          Reads V0s and tracks, writes out cascade vertices
+//    Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
+//------------------------------------------------------------------
+
+#include "TObject.h"
+
+class AliESDEvent;
+class AliESDv0;
+class AliExternalTrackParam;
+
+//_____________________________________________________________________________
+class AliLightCascadeVertexer : public TObject {
+public:
+  AliLightCascadeVertexer();
+  void SetCuts(const Double_t cuts[8]);
+  static void SetDefaultCuts(const Double_t cuts[8]);
+
+  Int_t V0sTracks2CascadeVertices(AliESDEvent *event);
+  Double_t Det(Double_t a00, Double_t a01, Double_t a10, Double_t a11) const;
+  Double_t Det(Double_t a00,Double_t a01,Double_t a02,
+              Double_t a10,Double_t a11,Double_t a12,
+              Double_t a20,Double_t a21,Double_t a22) const;
+
+  Double_t PropagateToDCA(AliESDv0 *vtx,AliExternalTrackParam *trk,Double_t b);
+
+  void GetCuts(Double_t cuts[8]) const;
+  static void GetDefaultCuts(Double_t cuts[8]);
+    
+    static void SetDefaultMaxEta(Double_t lMaxEta);
+    static void SetDefaultMinClusters(Double_t lMaxEta);
+    void SetMaxEta(Double_t lMaxEta);
+    void SetMinClusters(Double_t lMaxEta);
+
+private:
+  static
+  Double_t fgChi2max;   // maximal allowed chi2 
+  static
+  Double_t fgDV0min;    // min. allowed V0 impact parameter
+  static
+  Double_t fgMassWin;   // window around the Lambda mass
+  static
+  Double_t fgDBachMin;  // min. allowed bachelor impact parameter
+  static
+  Double_t fgDCAmax;    // maximal allowed DCA between the V0 and the track 
+  static
+  Double_t fgCPAmin;    // minimal allowed cosine of the cascade pointing angle
+  static
+  Double_t fgRmin, fgRmax;// max & min radii of the fiducial volume
+    static Double_t fgMaxEta;       // maximum eta value for track pre-selection
+    static Double_t fgMinClusters;  // minimum single-track clusters value (>=)
+  
+  Double_t fChi2max;    // maximal allowed chi2 
+  Double_t fDV0min;     // min. allowed V0 impact parameter
+  Double_t fMassWin;    // window around the Lambda mass
+  Double_t fDBachMin;   // min. allowed bachelor impact parameter
+  Double_t fDCAmax;     // maximal allowed DCA between the V0 and the track 
+  Double_t fCPAmin;     // minimal allowed cosine of the cascade pointing angle
+  Double_t fRmin, fRmax;// max & min radii of the fiducial volume
+    Double_t fMaxEta;       // maximum eta value for track pre-selection
+    Double_t fMinClusters;  // minimum single-track clusters value (>=)
+  
+  ClassDef(AliLightCascadeVertexer,3)  // cascade verterxer 
+};
+
+inline AliLightCascadeVertexer::AliLightCascadeVertexer() :
+  TObject(),
+  fChi2max(fgChi2max), 
+  fDV0min(fgDV0min),
+  fMassWin(fgMassWin),
+  fDBachMin(fgDBachMin),
+  fDCAmax(fgDCAmax),
+  fCPAmin(fgCPAmin), 
+  fRmin(fgRmin),
+fRmax(fgRmax),
+fMaxEta(fgMaxEta),
+fMinClusters(fgMinClusters)
+{
+}
+
+inline void AliLightCascadeVertexer::SetCuts(const Double_t cuts[8]) {
+  fChi2max=cuts[0]; 
+  fDV0min=cuts[1];   fMassWin=cuts[2]; fDBachMin=cuts[3];
+  fDCAmax=cuts[4];   fCPAmin=cuts[5];
+  fRmin=cuts[6];     fRmax=cuts[7]; 
+}
+
+inline void AliLightCascadeVertexer::SetDefaultCuts(const Double_t cuts[8]) {
+  fgChi2max=cuts[0]; 
+  fgDV0min=cuts[1];   fgMassWin=cuts[2]; fgDBachMin=cuts[3];
+  fgDCAmax=cuts[4];   fgCPAmin=cuts[5];
+  fgRmin=cuts[6];     fgRmax=cuts[7]; 
+}
+
+inline void AliLightCascadeVertexer::GetCuts(Double_t cuts[8]) const {
+  cuts[0]=fChi2max; 
+  cuts[1]=fDV0min;   cuts[2]=fMassWin;  cuts[3]=fDBachMin;
+  cuts[4]=fDCAmax;   cuts[5]=fCPAmin;
+  cuts[6]=fRmin;     cuts[7]=fRmax; 
+}
+
+inline void AliLightCascadeVertexer::GetDefaultCuts(Double_t cuts[8]) {
+  cuts[0]=fgChi2max; 
+  cuts[1]=fgDV0min;   cuts[2]=fgMassWin;  cuts[3]=fgDBachMin;
+  cuts[4]=fgDCAmax;   cuts[5]=fgCPAmin;
+  cuts[6]=fgRmin;     cuts[7]=fgRmax; 
+}
+
+inline void AliLightCascadeVertexer::SetDefaultMaxEta(Double_t lMaxEta) {
+    fgMaxEta = lMaxEta;
+}
+inline void AliLightCascadeVertexer::SetDefaultMinClusters(Double_t lMinClusters) {
+    fgMinClusters = lMinClusters;
+}
+inline void AliLightCascadeVertexer::SetMaxEta(Double_t lMaxEta) {
+    fMaxEta = lMaxEta;
+}
+inline void AliLightCascadeVertexer::SetMinClusters(Double_t lMinClusters) {
+    fMinClusters = lMinClusters;
+}
+
+#endif
+
diff --git a/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightV0vertexer.cxx b/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightV0vertexer.cxx
new file mode 100644 (file)
index 0000000..b3c2a4e
--- /dev/null
@@ -0,0 +1,170 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+//-------------------------------------------------------------------------
+//               Implementation of the V0 vertexer class
+//                  reads tracks writes out V0 vertices
+//                      fills the ESD with the V0s
+//     Origin: Iouri Belikov, IPHC, Strasbourg, Jouri.Belikov@cern.ch
+//-------------------------------------------------------------------------
+//                Adapted for more customizability and
+//                   lightweight operation in Pb-Pb
+//
+//          This is still being tested! Use at your own risk!
+//-------------------------------------------------------------------------
+
+#include "AliESDEvent.h"
+#include "AliESDv0.h"
+#include "AliLightV0vertexer.h"
+
+ClassImp(AliLightV0vertexer)
+
+
+//A set of very loose cuts
+Double_t AliLightV0vertexer::fgChi2max=33.; //max chi2
+Double_t AliLightV0vertexer::fgDNmin=0.05;  //min imp parameter for the 1st daughter
+Double_t AliLightV0vertexer::fgDPmin=0.05;  //min imp parameter for the 2nd daughter
+Double_t AliLightV0vertexer::fgDCAmax=1.5;  //max DCA between the daughter tracks
+Double_t AliLightV0vertexer::fgCPAmin=0.9;  //min cosine of V0's pointing angle
+Double_t AliLightV0vertexer::fgRmin=0.2;    //min radius of the fiducial volume
+Double_t AliLightV0vertexer::fgRmax=200.;   //max radius of the fiducial volume
+
+Double_t AliLightV0vertexer::fgMaxEta=0.8;        //max |eta|
+Double_t AliLightV0vertexer::fgMinClusters=70;   //min clusters (>=)
+
+Int_t AliLightV0vertexer::Tracks2V0vertices(AliESDEvent *event) {
+    //--------------------------------------------------------------------
+    //This function reconstructs V0 vertices
+    //--------------------------------------------------------------------
+    
+    const AliESDVertex *vtxT3D=event->GetPrimaryVertex();
+    
+    Double_t xPrimaryVertex=vtxT3D->GetX();
+    Double_t yPrimaryVertex=vtxT3D->GetY();
+    Double_t zPrimaryVertex=vtxT3D->GetZ();
+    
+    Int_t nentr=event->GetNumberOfTracks();
+    Double_t b=event->GetMagneticField();
+    
+    if (nentr<2) return 0;
+    
+    TArrayI neg(nentr);
+    TArrayI pos(nentr);
+    
+    Int_t nneg=0, npos=0, nvtx=0;
+    
+    Int_t i;
+    for (i=0; i<nentr; i++) {
+        AliESDtrack *esdTrack=event->GetTrack(i);
+        ULong_t status=esdTrack->GetStatus();
+        
+        //if ((status&AliESDtrack::kITSrefit)==0)//not to accept the ITS SA tracks
+        if ((status&AliESDtrack::kTPCrefit)==0) continue;
+        
+        //Track pre-selection: clusters
+        if (esdTrack->GetTPCNcls() < fMinClusters ) continue;
+        
+        Double_t d=esdTrack->GetD(xPrimaryVertex,yPrimaryVertex,b);
+        if (TMath::Abs(d)<fDPmin) continue;
+        if (TMath::Abs(d)>fRmax) continue;
+        
+        if (esdTrack->GetSign() < 0.) neg[nneg++]=i;
+        else pos[npos++]=i;
+    }
+    
+    
+    for (i=0; i<nneg; i++) {
+        Int_t nidx=neg[i];
+        AliESDtrack *ntrk=event->GetTrack(nidx);
+        
+        for (Int_t k=0; k<npos; k++) {
+            Int_t pidx=pos[k];
+            AliESDtrack *ptrk=event->GetTrack(pidx);
+            
+            //Track pre-selection: clusters
+            if (ptrk->GetTPCNcls() < fMinClusters ) continue;
+            
+            if (TMath::Abs(ntrk->GetD(xPrimaryVertex,yPrimaryVertex,b))<fDNmin)
+                if (TMath::Abs(ptrk->GetD(xPrimaryVertex,yPrimaryVertex,b))<fDNmin) continue;
+            
+            Double_t xn, xp, dca=ntrk->GetDCA(ptrk,b,xn,xp);
+            if (dca > fDCAmax) continue;
+            if ((xn+xp) > 2*fRmax) continue;
+            if ((xn+xp) < 2*fRmin) continue;
+            
+            AliExternalTrackParam nt(*ntrk), pt(*ptrk);
+            Bool_t corrected=kFALSE;
+            if ((nt.GetX() > 3.) && (xn < 3.)) {
+                //correct for the beam pipe material
+                corrected=kTRUE;
+            }
+            if ((pt.GetX() > 3.) && (xp < 3.)) {
+                //correct for the beam pipe material
+                corrected=kTRUE;
+            }
+            if (corrected) {
+                dca=nt.GetDCA(&pt,b,xn,xp);
+                if (dca > fDCAmax) continue;
+                if ((xn+xp) > 2*fRmax) continue;
+                if ((xn+xp) < 2*fRmin) continue;
+            }
+            
+            nt.PropagateTo(xn,b); pt.PropagateTo(xp,b);
+            
+            //select maximum eta range (after propagation)
+            if (TMath::Abs(nt.Eta())>fMaxEta) continue;
+            if (TMath::Abs(pt.Eta())>fMaxEta) continue;
+            
+            AliESDv0 vertex(nt,nidx,pt,pidx);
+            if (vertex.GetChi2V0() > fChi2max) continue;
+            
+            Double_t x=vertex.Xv(), y=vertex.Yv();
+            Double_t r2=x*x + y*y;
+            if (r2 < fRmin*fRmin) continue;
+            if (r2 > fRmax*fRmax) continue;
+            
+            Float_t cpa=vertex.GetV0CosineOfPointingAngle(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex);
+            
+            //Simple cosine cut (no pt dependence for now)
+            if (cpa < fCPAmin) continue;
+            
+            vertex.SetDcaV0Daughters(dca);
+            vertex.SetV0CosineOfPointingAngle(cpa);
+            vertex.ChangeMassHypothesis(kK0Short);
+            
+            event->AddV0(&vertex);
+            
+            nvtx++;
+        }
+    }
+    
+    Info("Tracks2V0vertices","Number of reconstructed V0 vertices: %d",nvtx);
+    
+    return nvtx;
+}
+
+
+
+
+
+
+
+
+
+
+
+
+
+
diff --git a/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightV0vertexer.h b/PWGLF/STRANGENESS/Cascades/lightvertexers/AliLightV0vertexer.h
new file mode 100644 (file)
index 0000000..b23e277
--- /dev/null
@@ -0,0 +1,122 @@
+#ifndef AliLightV0vertexer_H
+#define AliLightV0vertexer_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//------------------------------------------------------------------
+//                    V0 Vertexer Class
+//            reads tracks writes out V0 vertices
+//   Origin: Iouri Belikov, IReS, Strasbourg, Jouri.Belikov@cern.ch
+//------------------------------------------------------------------
+
+#include "TObject.h"
+
+class TTree;
+class AliESDEvent;
+
+//_____________________________________________________________________________
+class AliLightV0vertexer : public TObject {
+public:
+    AliLightV0vertexer();
+    void SetCuts(const Double_t cuts[7]);
+    static void SetDefaultCuts(const Double_t cuts[7]);
+    
+    Int_t Tracks2V0vertices(AliESDEvent *event);
+    
+    void GetCuts(Double_t cuts[7]) const;
+    static void GetDefaultCuts(Double_t cuts[7]);
+    
+    static void SetDefaultMaxEta(Double_t lMaxEta);
+    static void SetDefaultMinClusters(Double_t lMaxEta);
+    void SetMaxEta(Double_t lMaxEta);
+    void SetMinClusters(Double_t lMaxEta);
+
+    
+private:
+    static
+    Double_t fgChi2max;      // maximal allowed chi2
+    static
+    Double_t fgDNmin;        // min allowed impact parameter for the 1st daughter
+    static
+    Double_t fgDPmin;        // min allowed impact parameter for the 2nd daughter
+    static
+    Double_t fgDCAmax;       // maximal allowed DCA between the daughter tracks
+    static
+    Double_t fgCPAmin;       // minimal allowed cosine of V0's pointing angle
+    static
+    Double_t fgRmin, fgRmax; // max & min radii of the fiducial volume
+    
+    static Double_t fgMaxEta;       // maximum eta value for track pre-selection
+    static Double_t fgMinClusters;  // minimum single-track clusters value (>=)
+    
+    Double_t fChi2max;      // maximal allowed chi2
+    Double_t fDNmin;        // min allowed impact parameter for the 1st daughter
+    Double_t fDPmin;        // min allowed impact parameter for the 2nd daughter
+    Double_t fDCAmax;       // maximal allowed DCA between the daughter tracks
+    Double_t fCPAmin;       // minimal allowed cosine of V0's pointing angle
+    Double_t fRmin, fRmax;  // max & min radii of the fiducial volume
+    
+    Double_t fMaxEta;       // maximum eta value for track pre-selection
+    Double_t fMinClusters;  // minimum single-track clusters value (>=)
+    
+    ClassDef(AliLightV0vertexer,3)  // V0 verterxer
+};
+
+inline AliLightV0vertexer::AliLightV0vertexer() :
+TObject(),
+fChi2max(fgChi2max),
+fDNmin(fgDNmin),
+fDPmin(fgDPmin),
+fDCAmax(fgDCAmax),
+fCPAmin(fgCPAmin),
+fRmin(fgRmin),
+fRmax(fgRmax),
+fMaxEta(fgMaxEta),
+fMinClusters(fgMinClusters)
+{
+}
+
+inline void AliLightV0vertexer::SetCuts(const Double_t cuts[7]) {
+    fChi2max=cuts[0];
+    fDNmin=cuts[1];   fDPmin=cuts[2];
+    fDCAmax=cuts[3];  fCPAmin=cuts[4];
+    fRmin=cuts[5];    fRmax=cuts[6];
+}
+
+inline void AliLightV0vertexer::SetDefaultCuts(const Double_t cuts[7]) {
+    fgChi2max=cuts[0];
+    fgDNmin=cuts[1];   fgDPmin=cuts[2];
+    fgDCAmax=cuts[3];  fgCPAmin=cuts[4];
+    fgRmin=cuts[5];    fgRmax=cuts[6];
+}
+
+inline void AliLightV0vertexer::GetCuts(Double_t cuts[7]) const {
+    cuts[0]=fChi2max;
+    cuts[1]=fDNmin;   cuts[2]=fDPmin;
+    cuts[3]=fDCAmax;  cuts[4]=fCPAmin;
+    cuts[5]=fRmin;    cuts[6]=fRmax;
+}
+
+inline void AliLightV0vertexer::GetDefaultCuts(Double_t cuts[7]) {
+    cuts[0]=fgChi2max;
+    cuts[1]=fgDNmin;   cuts[2]=fgDPmin;
+    cuts[3]=fgDCAmax;  cuts[4]=fgCPAmin;
+    cuts[5]=fgRmin;    cuts[6]=fgRmax; 
+}
+
+inline void AliLightV0vertexer::SetDefaultMaxEta(Double_t lMaxEta) {
+    fgMaxEta = lMaxEta;
+}
+inline void AliLightV0vertexer::SetDefaultMinClusters(Double_t lMinClusters) {
+    fgMinClusters = lMinClusters;
+}
+inline void AliLightV0vertexer::SetMaxEta(Double_t lMaxEta) {
+    fMaxEta = lMaxEta;
+}
+inline void AliLightV0vertexer::SetMinClusters(Double_t lMinClusters) {
+    fMinClusters = lMinClusters;
+}
+
+#endif
+
+