]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
AliVertexerTracks used to estimate the position of the primary vertex (Yu.Belikov)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 16 May 2006 20:52:56 +0000 (20:52 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 16 May 2006 20:52:56 +0000 (20:52 +0000)
STEER/AliESD.cxx
STEER/AliESD.h
STEER/AliReconstruction.cxx
STEER/AliVertexerTracks.cxx
STEER/AliVertexerTracks.h

index 16143b8a1019c82c592b65d24996f39836d9b5b1..60beb8b8ae0125b1da7bfc660f00254118cb392c 100644 (file)
@@ -42,8 +42,9 @@ AliESD::AliESD():
   fZDCEMEnergy(0),
   fZDCParticipants(0),
   fT0zVertex(0),
+  fSPDVertex(0),
+  fPrimaryVertex(0),
   fT0timeStart(0),
-  fPrimaryVertex(),
   fTracks("AliESDtrack",15000),
   fHLTConfMapTracks("AliESDHLTtrack",25000),
   fHLTHoughTracks("AliESDHLTtrack",15000),
@@ -73,6 +74,8 @@ AliESD::~AliESD()
   //
   // Standard destructor
   //
+  delete fSPDVertex;
+  delete fPrimaryVertex;
   fTracks.Delete();
   fHLTConfMapTracks.Delete();
   fHLTHoughTracks.Delete();
@@ -124,7 +127,8 @@ void AliESD::Reset()
   fZDCParticipants=0;
   fT0zVertex=0;
   fT0timeStart = 0;
-  fPrimaryVertex.Reset();
+  delete fSPDVertex; fSPDVertex=0;
+  delete fPrimaryVertex; fPrimaryVertex=0;
   fTracks.Clear();
   fHLTConfMapTracks.Clear();
   fHLTHoughTracks.Clear();
@@ -154,9 +158,9 @@ void AliESD::Print(Option_t *) const
         GetTriggerMask(),
         GetMagneticField() );
   printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
-        fPrimaryVertex.GetXv(), fPrimaryVertex.GetXRes(),
-        fPrimaryVertex.GetYv(), fPrimaryVertex.GetYRes(),
-        fPrimaryVertex.GetZv(), fPrimaryVertex.GetZRes());
+        fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
+        fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
+        fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
   printf("Event from reconstruction version %d \n",fRecoVersion);
   printf("Number of tracks: \n");
   printf("                 charged   %d\n", GetNumberOfTracks());
index 1ac82984bf8453f565d25ebaec51e10588ef0f52..44fccd1bf7e862a8177403b48a5a5c571da14a46 100644 (file)
@@ -131,9 +131,14 @@ public:
   }
     
   void SetVertex(const AliESDVertex* vertex) {
-    new(&fPrimaryVertex) AliESDVertex(*vertex);
+     fSPDVertex=new AliESDVertex(*vertex);
   }
-  const AliESDVertex* GetVertex() const {return &fPrimaryVertex;};
+  const AliESDVertex* GetVertex() const {return fSPDVertex;};
+
+  void SetPrimaryVertex(const AliESDVertex* vertex) {
+     fPrimaryVertex=new AliESDVertex(*vertex);
+  }
+  const AliESDVertex* GetPrimaryVertex() const {return fPrimaryVertex;};
 
   Int_t  GetEventNumber() const {return fEventNumber;}
   Int_t  GetRunNumber() const {return fRunNumber;}
@@ -216,10 +221,12 @@ protected:
   Int_t        fZDCParticipants; // number of participants estimated by the ZDC
 
   Float_t      fT0zVertex;       // vertex z position estimated by the START
+  AliESDVertex *fSPDVertex;      // Primary vertex estimated by the SPD
+  AliESDVertex *fPrimaryVertex;  // Primary vertex estimated using ESD tracks
+
   Float_t      fT0timeStart;     // interaction time estimated by the START
   Float_t      fT0time[24];      // best TOF on each START PMT
   Float_t      fT0amplitude[24]; // number of particles(MIPs) on each START PMT
-  AliESDVertex fPrimaryVertex;   // Primary vertex estimated by the ITS
 
   TClonesArray fTracks;          // ESD tracks
   TClonesArray fHLTConfMapTracks;// HLT ESD tracks from Conformal Mapper method
@@ -240,7 +247,7 @@ protected:
  
   AliESDFMD *  fESDFMD; // FMD object containing rough multiplicity
 
-  ClassDef(AliESD,10)  //ESD class 
+  ClassDef(AliESD,11)  //ESD class 
 };
 #endif 
 
index a49ca602e10c5d2977a878cccc77828bddd86acf..ae1144b104d12a96ef9769a22e1f82cb6a05741c 100644 (file)
 #include "AliESDVertex.h"
 #include "AliTracker.h"
 #include "AliVertexer.h"
+#include "AliVertexerTracks.h"
 #include "AliHeader.h"
 #include "AliGenEventHeader.h"
 #include "AliPID.h"
@@ -585,6 +586,8 @@ Bool_t AliReconstruction::Run(const char* input,
 
   gROOT->cd();
 
+  AliVertexerTracks tVertexer;
+
   // loop over events
   if (fRawReader) fRawReader->RewindEvents();
   
@@ -683,6 +686,8 @@ Bool_t AliReconstruction::Run(const char* input,
       }
     }
 
+    esd->SetPrimaryVertex(tVertexer.FindVertex(esd));
+
     // write ESD
     tree->Fill();
     // write HLT ESD
index c3744abb3b02ab3ecfa12ac8b6fc090b6949f3cb..12bbf299ef516d84630a4cfdcf700828e3072915 100644 (file)
@@ -30,6 +30,7 @@
 //---- AliRoot headers -----
 #include "AliStrLine.h"
 #include "AliVertexerTracks.h"
+#include "AliESD.h"
 #include "AliESDtrack.h"
 
 ClassImp(AliVertexerTracks)
@@ -116,6 +117,43 @@ AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) {
   if(fAlgo==5)  VertexFinder(0);
   return &fVert;
 }
+
+//----------------------------------------------------------------------------
+AliESDVertex* AliVertexerTracks::FindVertex(const AliESD *event) {
+//
+// This is a simple wrapping (by Jouri.Belikov@cern.ch) over the original
+// code by the authors of this class.
+//
+  Int_t nt=event->GetNumberOfTracks(), nacc=0;
+  while (nt--) {
+    AliESDtrack *t=event->GetTrack(nt);
+    if ((t->GetStatus()&AliESDtrack::kITSrefit)==0) continue;
+    fTrkArray.AddLast(t);
+    nacc++;
+  }
+
+  // get tracks and propagate them to initial vertex position
+  if(nacc < fMinTracks) {
+    printf("TooFewTracks\n");
+    Double_t vtx[3]={0,0,0};
+    fVert.SetXYZ(vtx);
+    fVert.SetDispersion(999);
+    fVert.SetNContributors(-5);
+  } else 
+    switch (fAlgo) {
+       case 1: StrLinVertexFinderMinDist(1); break;
+       case 2: StrLinVertexFinderMinDist(0); break;
+       case 3: HelixVertexFinder();          break;
+       case 4: VertexFinder(1);              break;
+       case 5: VertexFinder(0);              break;
+       default: printf("Wrong algorithm\n"); break;  
+    }
+
+  fTrkArray.Clear();
+  return &fVert;
+}
+
+
 //----------------------------------------------------------------------------
 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) {
 //
index 8f3bc6b782b6d39168e96c511dcac13eecf9c6cd..4f12fb11046890cc7d6ca811d6133a2a7f831213 100644 (file)
  *                                                                           *
  *****************************************************************************/
 
-#include "AliVertex.h"
+#include "AliESDVertex.h"
 #include "AliTracker.h"
 
 #include <TObjArray.h>
 
 class TTree; 
-class AliVertex; 
 class AliESD;
 
 class AliVertexerTracks : public TObject {
@@ -37,6 +36,8 @@ class AliVertexerTracks : public TObject {
   AliVertexerTracks(Double_t xStart, Double_t yStart); 
   virtual ~AliVertexerTracks();
 
+  AliESDVertex *FindVertex(const AliESD *event);
+
   // computes the vertex from the set of tracks in the tree
   AliVertex* VertexForSelectedTracks(TTree *trkTree);
   AliVertex* VertexForSelectedTracks(TObjArray *trkArray);
@@ -61,7 +62,7 @@ class AliVertexerTracks : public TObject {
   static Double_t GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0);
   static Double_t GetDeterminant3X3(Double_t matr[][3]);
 
-  AliVertex fVert;            // vertex after vertex finder
+  AliESDVertex fVert;         // vertex after vertex finder
   Double_t  fNominalPos[2];   // initial knowledge on vertex position
   Int_t     fMinTracks;       // minimum number of tracks
   TObjArray fTrkArray;        // array with tracks to be processed
@@ -81,7 +82,7 @@ class AliVertexerTracks : public TObject {
   //         approximated as straight lines 
 
 
-  ClassDef(AliVertexerTracks,1) // 3D Vertexing with ESD tracks 
+  ClassDef(AliVertexerTracks,2) // 3D Vertexing with ESD tracks 
 };
 
 #endif