]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Introducing the ESD friend classes for storing complementary to ESD information in...
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 23 Apr 2006 05:50:15 +0000 (05:50 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 23 Apr 2006 05:50:15 +0000 (05:50 +0000)
STEER/AliESDfriend.cxx [new file with mode: 0644]
STEER/AliESDfriend.h [new file with mode: 0644]
STEER/AliESDfriendTrack.cxx [new file with mode: 0644]
STEER/AliESDfriendTrack.h [new file with mode: 0644]
STEER/AliExternalTrackParam.cxx
STEER/AliExternalTrackParam.h
STEER/AliReconstruction.cxx
STEER/AliReconstruction.h
STEER/ESDLinkDef.h
STEER/libESD.pkg

diff --git a/STEER/AliESDfriend.cxx b/STEER/AliESDfriend.cxx
new file mode 100644 (file)
index 0000000..f8b9efa
--- /dev/null
@@ -0,0 +1,61 @@
+/**************************************************************************
+ * 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 AliESDfriend class
+//  This class contains some additional to the ESD information like
+//  the clusters associated to tracks.
+//      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
+//-------------------------------------------------------------------------
+
+#include "AliESDfriend.h"
+#include "AliESDfriendTrack.h"
+#include "AliESD.h"
+
+ClassImp(AliESDfriend)
+
+AliESDfriend::AliESDfriend(): TObject(), fTracks("AliESDfriendTrack",15000)
+{
+ //
+ // Default constructor
+ //
+}
+
+AliESDfriend::AliESDfriend(const AliESDfriend &f):TObject(f),fTracks(f.fTracks)
+{
+ //
+ // Copy constructor
+ //
+}
+
+AliESDfriend::AliESDfriend(const AliESD &event): TObject(event),
+fTracks("AliESDfriendTrack",event.GetNumberOfTracks()) {
+  //
+  // Extracts the additional info from the ESD
+  //
+  Int_t ntrk=event.GetNumberOfTracks();
+
+  for (Int_t i=0; i<ntrk; i++) {
+    const AliESDtrack *t=event.GetTrack(i);
+    new (fTracks[fTracks.GetEntriesFast()]) AliESDfriendTrack(*t); 
+  }
+}
+
+AliESDfriend::~AliESDfriend() {
+  //
+  // Destructor
+  //
+  fTracks.Delete();
+}
diff --git a/STEER/AliESDfriend.h b/STEER/AliESDfriend.h
new file mode 100644 (file)
index 0000000..b59da6a
--- /dev/null
@@ -0,0 +1,36 @@
+#ifndef ALIESDFRIEND_H
+#define ALIESDFRIEND_H
+
+//-------------------------------------------------------------------------
+//                     Class AliESDfriend
+//               This class contains ESD additions
+//       Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch 
+//-------------------------------------------------------------------------
+
+#include <TObject.h>
+#include <TClonesArray.h>
+
+class AliESD;
+class AliESDfriendTrack;
+
+//_____________________________________________________________________________
+class AliESDfriend : public TObject {
+public:
+  AliESDfriend();
+  AliESDfriend(const AliESDfriend &);
+  AliESDfriend(const AliESD &);
+  virtual ~AliESDfriend();
+
+  Int_t GetNumberOfTracks() const {return fTracks.GetEntriesFast();}
+  AliESDfriendTrack *GetTrack(Int_t i) const {
+    return (AliESDfriendTrack *)fTracks.UncheckedAt(i);
+  }
+
+protected:
+  TClonesArray fTracks;    // ESD friend tracks
+  ClassDef(AliESDfriend,1) // ESD friend
+};
+
+#endif
+
+
diff --git a/STEER/AliESDfriendTrack.cxx b/STEER/AliESDfriendTrack.cxx
new file mode 100644 (file)
index 0000000..72e7b79
--- /dev/null
@@ -0,0 +1,67 @@
+/**************************************************************************
+ * 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 AliESDfriendTrack class
+//  This class keeps complementary to the AliESDtrack information 
+//      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
+//-------------------------------------------------------------------------
+#include "AliTrackPointArray.h"
+#include "AliESDfriendTrack.h"
+#include "AliESD.h"
+
+ClassImp(AliESDfriendTrack)
+
+  AliESDfriendTrack::AliESDfriendTrack(): TObject(), f1P(0), fPoints(0) {
+  //
+  // Default constructor
+  //
+  Int_t i;
+  for (i=0; i<AliESDtrack::kMaxITScluster; i++) fITSindex[i]=-2;
+  for (i=0; i<AliESDtrack::kMaxTPCcluster; i++) fTPCindex[i]=-2;
+  for (i=0; i<AliESDtrack::kMaxTRDcluster; i++) fTRDindex[i]=-2;
+}
+
+AliESDfriendTrack::AliESDfriendTrack(const AliESDfriendTrack &t): 
+TObject(t),
+f1P(t.f1P),
+fPoints(0)
+{
+  //
+  // Copy constructor
+  //
+  Int_t i;
+  for (i=0; i<AliESDtrack::kMaxITScluster; i++) fITSindex[i]=t.fITSindex[i];
+  for (i=0; i<AliESDtrack::kMaxTPCcluster; i++) fTPCindex[i]=t.fTPCindex[i];
+  for (i=0; i<AliESDtrack::kMaxTRDcluster; i++) fTRDindex[i]=t.fTRDindex[i];
+  if (t.fPoints) fPoints=new AliTrackPointArray(*t.fPoints);
+}
+
+AliESDfriendTrack::AliESDfriendTrack(const AliESDtrack &t): 
+TObject(t),
+f1P(t.Get1P()),
+fPoints(0) 
+{
+  //
+  // Extracts the complementary info from the ESD track
+  //
+  t.GetITSclusters(fITSindex); 
+  t.GetTPCclusters(fTPCindex); 
+  t.GetTRDclusters(fTRDindex); 
+  const AliTrackPointArray *points=t.GetTrackPointArray();
+  if (points) fPoints=new AliTrackPointArray(*points);
+}
+
+AliESDfriendTrack::~AliESDfriendTrack() {delete fPoints;}
diff --git a/STEER/AliESDfriendTrack.h b/STEER/AliESDfriendTrack.h
new file mode 100644 (file)
index 0000000..69ee70f
--- /dev/null
@@ -0,0 +1,41 @@
+#ifndef ALIESDFRIENDTRACK_H
+#define ALIESDFRIENDTRACK_H
+
+//-------------------------------------------------------------------------
+//                     Class AliESDfriendTrack
+//               This class contains ESD track additions
+//       Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch 
+//-------------------------------------------------------------------------
+
+#include <TObject.h>
+#include "AliESDtrack.h"
+
+class AliTrackPointArray;
+
+//_____________________________________________________________________________
+class AliESDfriendTrack : public TObject {
+public:
+  AliESDfriendTrack();
+  AliESDfriendTrack(const AliESDfriendTrack &);
+  AliESDfriendTrack(const AliESDtrack &);
+  virtual ~AliESDfriendTrack();
+
+  Float_t Get1P() const {return f1P;}
+  const Int_t *GetITSindices() const {return fITSindex;}
+  const Int_t *GetTPCindices() const {return fTPCindex;}
+  const Int_t *GetTRDindices() const {return fTRDindex;}
+  const AliTrackPointArray *GetTrackPointArray() const {return fPoints;}
+
+protected:
+  Float_t f1P;                                  // 1/P (1/(GeV/c))
+  Int_t fITSindex[AliESDtrack::kMaxITScluster]; // indices of the ITS clusters 
+  Int_t fTPCindex[AliESDtrack::kMaxTPCcluster]; // indices of the TPC clusters
+  Int_t fTRDindex[AliESDtrack::kMaxTRDcluster]; // indices of the TRD clusters
+
+  AliTrackPointArray *fPoints; // Array of track space points in the global frame
+  ClassDef(AliESDfriendTrack,1) //ESD friend track
+};
+
+#endif
+
+
index 4db4b13877f801110b272a384bc5ffeeb91a4ce8..15e5ced21dcc30b986766e28202cd01cd00d0f6f 100644 (file)
@@ -91,10 +91,17 @@ Double_t AliExternalTrackParam::GetP() const {
   // This function returns the track momentum
   // Results for (nearly) straight tracks are meaningless !
   //---------------------------------------------------------------------
-  if (TMath::Abs(fP[4])<=0) return 0;
+  if (TMath::Abs(fP[4])<=0) return 1e+33;
   return TMath::Sqrt(1.+ fP[3]*fP[3])/TMath::Abs(fP[4]);
 }
 
+Double_t AliExternalTrackParam::Get1P() const {
+  //---------------------------------------------------------------------
+  // This function returns the 1/(track momentum)
+  //---------------------------------------------------------------------
+  return TMath::Abs(fP[4])/TMath::Sqrt(1.+ fP[3]*fP[3]);
+}
+
 //_______________________________________________________________________
 Double_t AliExternalTrackParam::GetD(Double_t x,Double_t y,Double_t b) const {
   //------------------------------------------------------------------
index b88fbcdf9a19bca9b6b968b05c801afb809e9241..599099704106d33fa65bd8cec1785295a6cb8158 100644 (file)
@@ -49,6 +49,7 @@ class AliExternalTrackParam: public TObject {
   Double_t GetAlpha() const {return fAlpha;}
   Double_t GetSign() const {return (fP[4]>0) ? 1 : -1;}
   Double_t GetP() const;
+  Double_t Get1P() const;
   Double_t GetD(Double_t xv, Double_t yv, Double_t b) const; 
   Double_t GetLinearD(Double_t xv, Double_t yv) const; 
   Bool_t CorrectForMaterial(Double_t d, Double_t x0, Double_t mass);
index 790090155215f314cb499dfc7dd5eb207438f57f..edf38ffa18177ce1f62161063c374f785e500dd0 100644 (file)
@@ -73,7 +73,7 @@
 //                                                                           //
 // Uniform/nonuniform field tracking switches (default: uniform field)       //
 //                                                                           //
-//   rec.SetUniformFieldTracking();  ( rec.SetNonuniformFieldTracking(); )   //
+//   rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
 //                                                                           //
 // The filling of additional ESD information can be steered by               //
 //                                                                           //
 #include "AliRawReaderDate.h"
 #include "AliRawReaderRoot.h"
 #include "AliESD.h"
+#include "AliESDfriend.h"
 #include "AliESDVertex.h"
 #include "AliTracker.h"
 #include "AliVertexer.h"
 #include "AliESDtrack.h"
 
 #include "AliRunTag.h"
-//#include "AliLHCTag.h"
 #include "AliDetectorTag.h"
 #include "AliEventTag.h"
 
@@ -153,17 +153,20 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename, const char* cdb
                                     const char* name, const char* title) :
   TNamed(name, title),
 
-  fRunLocalReconstruction("ALL"),
   fUniformField(kTRUE),
   fRunVertexFinder(kTRUE),
   fRunHLTTracking(kFALSE),
+  fStopOnError(kFALSE),
+  fWriteAlignmentData(kFALSE),
+  fWriteESDfriend(kFALSE),
+
+  fRunLocalReconstruction("ALL"),
   fRunTracking("ALL"),
   fFillESD("ALL"),
   fGAliceFileName(gAliceFilename),
   fInput(""),
   fFirstEvent(0),
   fLastEvent(-1),
-  fStopOnError(kFALSE),
   fCheckPointLevel(0),
   fOptions(),
   fLoadAlignFromCDB(kTRUE),
@@ -175,7 +178,6 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename, const char* cdb
   fVertexer(NULL),
 
   fAlignObjArray(NULL),
-  fWriteAlignmentData(kFALSE),
   fCDBUri(cdbUri)
 {
 // create reconstruction object with default parameters
@@ -192,17 +194,20 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename, const char* cdb
 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   TNamed(rec),
 
-  fRunLocalReconstruction(rec.fRunLocalReconstruction),
   fUniformField(rec.fUniformField),
   fRunVertexFinder(rec.fRunVertexFinder),
   fRunHLTTracking(rec.fRunHLTTracking),
+  fStopOnError(rec.fStopOnError),
+  fWriteAlignmentData(rec.fWriteAlignmentData),
+  fWriteESDfriend(rec.fWriteESDfriend),
+
+  fRunLocalReconstruction(rec.fRunLocalReconstruction),
   fRunTracking(rec.fRunTracking),
   fFillESD(rec.fFillESD),
   fGAliceFileName(rec.fGAliceFileName),
   fInput(rec.fInput),
   fFirstEvent(rec.fFirstEvent),
   fLastEvent(rec.fLastEvent),
-  fStopOnError(rec.fStopOnError),
   fCheckPointLevel(0),
   fOptions(),
   fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
@@ -214,7 +219,6 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fVertexer(NULL),
 
   fAlignObjArray(rec.fAlignObjArray),
-  fWriteAlignmentData(rec.fWriteAlignmentData),
   fCDBUri(rec.fCDBUri)
 {
 // copy constructor
@@ -562,6 +566,18 @@ Bool_t AliReconstruction::Run(const char* input,
   hlttree->Branch("ESD", "AliESD", &hltesd);
   delete esd; delete hltesd;
   esd = NULL; hltesd = NULL;
+
+  // create the file and tree with ESD additions
+  TFile *filef=0; TTree *treef=0; AliESDfriend *esdf=0;
+  if (fWriteESDfriend) {
+     filef = TFile::Open("AliESDfriends.root", "RECREATE");
+     if (!filef->IsOpen()) {
+        AliError("opening AliESDfriends.root failed");
+     }
+     treef = new TTree("esdFriendTree", "Tree with ESD friends");
+     treef->Branch("ESDfriend", "AliESDfriend", &esdf);
+  }
+
   gROOT->cd();
 
   // loop over events
@@ -658,6 +674,12 @@ Bool_t AliReconstruction::Run(const char* input,
     // write HLT ESD
     hlttree->Fill();
 
+    // write ESD friend
+    if (fWriteESDfriend) {
+       esdf=new AliESDfriend(*esd);
+       treef->Fill();  
+    }
+
     if (fCheckPointLevel > 0)  WriteESD(esd, "final"); 
  
     delete esd; delete hltesd;
@@ -671,6 +693,11 @@ Bool_t AliReconstruction::Run(const char* input,
   tree->Write();
   hlttree->Write();
 
+  if (fWriteESDfriend) {
+     filef->cd();
+     treef->Write(); delete treef; filef->Close(); delete filef; 
+  }
+
   // Create tags for the events in the ESD tree (the ESD tree is always present)
   // In case of empty events the tags will contain dummy values
   CreateTag(file);
index 28f217c0d46563afcabe92c20b46277535dbf448..fa4bdc4a758acddf5a04b19a342f60ab9987e157 100644 (file)
@@ -54,7 +54,6 @@ public:
 
   void           SetRunLocalReconstruction(const char* detectors) {
     fRunLocalReconstruction = detectors;};
-  void           SetRunVertexFinder(Bool_t run) {fRunVertexFinder = run;};
   void           SetRunTracking(const char* detectors) {
     fRunTracking = detectors;};
   void           SetFillESD(const char* detectors) {fFillESD = detectors;};
@@ -66,14 +65,18 @@ public:
   void           SetLoadAlignData(const char* detectors) 
     {fLoadAlignData = detectors;};
 
-   void SetUniformFieldTracking(){fUniformField=kTRUE;} 
-   void SetNonuniformFieldTracking(){fUniformField=kFALSE;} 
 
-  void           SetStopOnError(Bool_t stopOnError) 
-    {fStopOnError = stopOnError;}
+  //*** Global reconstruction flag setters
+  void SetUniformFieldTracking(Bool_t flag=kTRUE){fUniformField=flag;} 
+  void SetRunVertexFinder(Bool_t flag=kTRUE) {fRunVertexFinder=flag;};
+  void SetRunHLTTracking(Bool_t flag=kTRUE) {fRunHLTTracking=flag;};
+  void SetStopOnError(Bool_t flag=kTRUE) {fStopOnError=flag;}
+  void SetWriteAlignmentData(Bool_t flag=kTRUE){fWriteAlignmentData=flag;}
+  void SetWriteESDfriend(Bool_t flag=kTRUE){fWriteESDfriend=flag;}
+
+                  
   void           SetCheckPointLevel(Int_t checkPointLevel)
     {fCheckPointLevel = checkPointLevel;}
-                  
   // CDB storage activation
   void InitCDBStorage();
   void SetDefaultStorage(const char* uri);
@@ -96,7 +99,6 @@ public:
   Bool_t         Run(Int_t firstEvent, Int_t lastEvent = -1)
     {return Run(NULL, firstEvent, lastEvent);};
 
-  void SetWriteAlignmentData(){fWriteAlignmentData=kTRUE;}
 
 private:
   Bool_t         RunLocalReconstruction(const TString& detectors);
@@ -124,18 +126,22 @@ private:
   void           WriteAlignmentData(AliESD* esd);
 
 
-
-  TString        fRunLocalReconstruction; // run the local reconstruction for these detectors
+  //*** Global reconstruction flags *******************
   Bool_t         fUniformField;       // uniform field tracking flag
   Bool_t         fRunVertexFinder;    // run the vertex finder
   Bool_t         fRunHLTTracking;     // run the HLT tracking
+  Bool_t         fStopOnError;        // stop or continue on errors
+  Bool_t         fWriteAlignmentData; // write track space-points flag
+  Bool_t         fWriteESDfriend;     // write ESD friend flag
+
+
+  TString        fRunLocalReconstruction; // run the local reconstruction for these detectors
   TString        fRunTracking;        // run the tracking for these detectors
   TString        fFillESD;            // fill ESD for these detectors
   TString        fGAliceFileName;     // name of the galice file
   TString        fInput;              // name of input file or directory
   Int_t          fFirstEvent;         // index of first event to be reconstr.
   Int_t          fLastEvent;          // index of last event to be reconstr.
-  Bool_t         fStopOnError;        // stop or continue on errors
   Int_t          fCheckPointLevel;    // level of ESD check points
   TObjArray      fOptions;            // options for reconstructor objects
   Bool_t         fLoadAlignFromCDB;   // Load alignment data from CDB and apply it to geometry or not
@@ -152,11 +158,10 @@ private:
   AliTracker*    fTracker[fgkNDetectors];  //! trackers
 
   TObjArray*    fAlignObjArray;      // array with the alignment objects to be applied to the geometry
-  Bool_t         fWriteAlignmentData; // write track space-points flag
 
   TString       fCDBUri;             // Uri of the default CDB storage
 
-  ClassDef(AliReconstruction, 6)      // class for running the reconstruction
+  ClassDef(AliReconstruction, 7)      // class for running the reconstruction
 };
 
 #endif
index 7139db91aa8a585991f4ca0f0f3451fa062ac959..abf00245d64e1080dd9cce1d6609fa9ac0379b74 100644 (file)
@@ -11,7 +11,9 @@
 #pragma link C++ enum   AliLog::EType_t;
 
 #pragma link C++ class  AliESD+;
+#pragma link C++ class  AliESDfriend+;
 #pragma link C++ class  AliESDtrack+;
+#pragma link C++ class  AliESDfriendTrack+;
 #pragma link C++ class  AliESDMuonTrack+;
 #pragma link C++ class  AliESDPmdTrack+;
 #pragma link C++ class  AliESDTrdTrack+;
index 6ca25f0fe83782bda2c1b70f03f846a4bb08a838..c377e44599db7fdc083fe4956f37597c2c3f2d5e 100644 (file)
@@ -1,5 +1,5 @@
-SRCS = AliESD.cxx \
-       AliESDtrack.cxx \
+SRCS = AliESD.cxx AliESDfriend.cxx\
+       AliESDtrack.cxx AliESDfriendTrack.cxx\
        AliESDMuonTrack.cxx AliESDPmdTrack.cxx AliESDTrdTrack.cxx AliESDHLTtrack.cxx \
        AliESDv0.cxx AliESDcascade.cxx AliVertex.cxx AliESDVertex.cxx \
        AliESDpid.cxx AliESDkink.cxx AliESDV0MI.cxx \