From 549b0ea14d028d11e3e1e234daf4ac1c00fccd81 Mon Sep 17 00:00:00 2001
From: Rene Brun <Rene.Brun@cern.ch>
Date: Fri, 5 Oct 2001 16:49:40 +0000
Subject: [PATCH] Modify the test program to illustrate the new classes TRef
 and TRefArray.

The following members have been added to Event.h
   TRefArray     *fHighPt;       //array of High Pt tracks only
   TRefArray     *fMuons;        //array of Muon tracks only
   TRef           fLastTrack;    //pointer to last track

In MainEvent (write part) the two arrays are filled in the loop
making the tracks. fLastTrack is set directly in AddTrack.

Split mode can be used to store in separate branches the original
tracks and the two selected track arrays. Only one copy of a track
is stored in the file. When the file is read back, one can read
in any order the TClonesArray of tracks or the other arrays.
If the array fMuons is read and the TClonesArray branches are not
read, fMuons->At(i) will return 0. When the TClonesArray is read,
fMuons->at(i) will return a pointer to the selected track number i
in the TRefArray fMuons. This will point directly to one of the
tracks in the TClonesArray.

Note that the branches fHighPt, fMuons and fLastTrack could also
be stored in separate files (may be friends of the master file).


git-svn-id: http://root.cern.ch/svn/root/trunk@2987 27541ba8-7e3a-0410-8455-c3a389f83636
---
 test/Event.cxx     |  30 +++++++----
 test/Event.h       | 124 ++++++++++++++++++++++++---------------------
 test/MainEvent.cxx |  25 +++++++--
 3 files changed, 106 insertions(+), 73 deletions(-)

diff --git a/test/Event.cxx b/test/Event.cxx
index 3a8125b7458..68dd15c75e8 100644
--- a/test/Event.cxx
+++ b/test/Event.cxx
@@ -1,4 +1,4 @@
-// @(#)root/test:$Name:  $:$Id: Event.cxx,v 1.9 2000/11/21 21:06:53 brun Exp $
+// @(#)root/test:$Name:  $:$Id: Event.cxx,v 1.10 2000/11/23 10:22:10 brun Exp $
 // Author: Rene Brun   19/08/96
 
 ////////////////////////////////////////////////////////////////////////
@@ -14,11 +14,15 @@
 //        Int_t          fNvertex;
 //        UInt_t         fFlag;
 //        Float_t        fTemperature;
+//        Int_t          fMeasures[10];
+//        Float_t        fMatrix[4][4];
+//        Float_t       *fClosestDistance; //[fNvertex] indexed array! 
 //        EventHeader    fEvtHdr;
 //        TClonesArray  *fTracks;
+//        TRefArray     *fHighPt;            //array of High Pt tracks only
+//        TRefArray     *fMuons;             //array of Muon tracks only
+//        TRef           fLastTrack;         //pointer to last track
 //        TH1F          *fH;
-//        Float_t        fMatrix[4][4];
-//        Float_t       *fClosestDistance; //[fNvertex] indexed array! 
 //
 //   The EventHeader class has 3 data members (integers):
 //     public:
@@ -85,6 +89,8 @@ Event::Event()
 
    if (!fgTracks) fgTracks = new TClonesArray("Track", 1000);
    fTracks = fgTracks;
+   fHighPt = new TRefArray;
+   fMuons  = new TRefArray;
    fNtrack = 0;
    fH      = 0;
    Int_t i0,i1;
@@ -102,13 +108,14 @@ Event::~Event()
 {
    Clear();
    if (fH == fgHist) fgHist = 0;
-   delete fH;
-   fH = 0;
+   delete fH; fH = 0;
+   delete fHighPt; fHighPt = 0;
+   delete fMuons;  fMuons = 0;
    delete [] fClosestDistance;
 }
 
 //______________________________________________________________________________
-void Event::AddTrack(Float_t random)
+Track *Event::AddTrack(Float_t random)
 {
    // Add a new track to the list of tracks for this event.
    // To avoid calling the very time consuming operator new for each track,
@@ -117,13 +124,17 @@ void Event::AddTrack(Float_t random)
    // otherwise the previous Track[i] will be overwritten.
 
    TClonesArray &tracks = *fTracks;
-   new(tracks[fNtrack++]) Track(random);
+   Track *track = new(tracks[fNtrack++]) Track(random);
+   fLastTrack = track;
+   return track;
 }
 
 //______________________________________________________________________________
 void Event::Clear(Option_t *option)
 {
    fTracks->Clear(option);
+   fHighPt->Delete();
+   fMuons->Delete();
 }
 
 //______________________________________________________________________________
@@ -131,6 +142,7 @@ void Event::Reset(Option_t *option)
 {
 // Static function to reset all static objects for this event
 //   fgTracks->Delete(option);
+
    delete fgTracks; fgTracks = 0;
    fgHist   = 0;
 }
@@ -153,7 +165,7 @@ void Event::SetMeasure(UChar_t which, Int_t what) {
 //______________________________________________________________________________
 void Event::SetRandomVertex() {
    // This delete is to test the relocation of variable length array
-   delete fClosestDistance;
+   if (fClosestDistance) delete [] fClosestDistance;
    if (!fNvertex) {
       fClosestDistance = 0;
       return;
@@ -176,7 +188,7 @@ Track::Track(Float_t random) : TObject()
    fPy = py;
    fPz = TMath::Sqrt(px*px+py*py);
    fRandom = 1000*random;
-   if (fRandom < 10) fMass2 = 0.08;
+   if (fRandom < 10) fMass2 = 0.106;
    else if (fRandom < 100) fMass2 = 0.8;
    else if (fRandom < 500) fMass2 = 4.5;
    else if (fRandom < 900) fMass2 = 8.9;
diff --git a/test/Event.h b/test/Event.h
index 48e60ecadb0..fd431ecc10c 100644
--- a/test/Event.h
+++ b/test/Event.h
@@ -11,14 +11,64 @@
 
 #include "TObject.h"
 #include "TClonesArray.h"
+#include "TRefArray.h"
+#include "TRef.h"
 #include "TH1.h"
 #include "TMath.h"
 
 #include <iostream.h>
 
-
 class TDirectory;
 
+class Track : public TObject {
+
+private:
+   Float_t      fPx;           //X component of the momentum
+   Float_t      fPy;           //Y component of the momentum
+   Float_t      fPz;           //Z component of the momentum
+   Float_t      fRandom;       //A random track quantity
+   Float_t      fMass2;        //The mass square of this particle
+   Float_t      fBx;           //X intercept at the vertex
+   Float_t      fBy;           //Y intercept at the vertex
+   Float_t      fMeanCharge;   //Mean charge deposition of all hits of this track
+   Float_t      fXfirst;       //X coordinate of the first point
+   Float_t      fXlast;        //X coordinate of the last point
+   Float_t      fYfirst;       //Y coordinate of the first point
+   Float_t      fYlast;        //Y coordinate of the last point
+   Float_t      fZfirst;       //Z coordinate of the first point
+   Float_t      fZlast;        //Z coordinate of the last point
+   Float_t      fCharge;       //Charge of this track
+   Float_t      fVertex[3];    //Track vertex position
+   Int_t        fNpoint;       //Number of points for this track
+   Short_t      fValid;        //Validity criterion
+
+public:
+   Track() { }
+   Track(Float_t random);
+   virtual ~Track() { }
+   Float_t       GetPx() const { return fPx; }
+   Float_t       GetPy() const { return fPy; }
+   Float_t       GetPz() const { return fPz; }
+   Float_t       GetPt() const { return TMath::Sqrt(fPx*fPx + fPy*fPy); }
+   Float_t       GetRandom() const { return fRandom; }
+   Float_t       GetBx() const { return fBx; }
+   Float_t       GetBy() const { return fBy; }
+   Float_t       GetMass2() const { return fMass2; }
+   Float_t       GetMeanCharge() const { return fMeanCharge; }
+   Float_t       GetXfirst() const { return fXfirst; }
+   Float_t       GetXlast()  const { return fXlast; }
+   Float_t       GetYfirst() const { return fYfirst; }
+   Float_t       GetYlast()  const { return fYlast; }
+   Float_t       GetZfirst() const { return fZfirst; }
+   Float_t       GetZlast()  const { return fZlast; }
+   Float_t       GetCharge() const { return fCharge; }
+   Float_t       GetVertex(Int_t i=0) {return fVertex[i];}
+   Int_t         GetNpoint() const { return fNpoint; }
+   Short_t       GetValid()  const { return fValid; }
+   virtual void  SetValid(Int_t valid=1) { fValid = valid; }
+
+   ClassDef(Track,1)  //A track segment
+};
 
 class EventHeader {
 
@@ -42,19 +92,22 @@ public:
 class Event : public TObject {
 
 private:
-   char           fType[20];
-   Int_t          fNtrack;
-   Int_t          fNseg;
+   char           fType[20];          //event type
+   Int_t          fNtrack;            //Number of tracks
+   Int_t          fNseg;              //Number of track segments
    Int_t          fNvertex;
    UInt_t         fFlag;
    Float_t        fTemperature;
-   EventHeader    fEvtHdr;
-   TClonesArray  *fTracks;            //->
-   TH1F          *fH;                 //->
    Int_t          fMeasures[10];
    Float_t        fMatrix[4][4];
    Float_t       *fClosestDistance;   //[fNvertex] 
-
+   EventHeader    fEvtHdr;
+   TClonesArray  *fTracks;            //->array with all tracks
+   TRefArray     *fHighPt;            //array of High Pt tracks only
+   TRefArray     *fMuons;             //array of Muon tracks only
+   TRef           fLastTrack;         //reference pointer to last track
+   TH1F          *fH;                 //->
+      
    static TClonesArray *fgTracks;
    static TH1F         *fgHist;
 
@@ -71,7 +124,7 @@ public:
    void          SetTemperature(Float_t t) { fTemperature = t; }
    void          SetType(char *type) {strcpy(fType,type);}
    void          SetHeader(Int_t i, Int_t run, Int_t date, Float_t random);
-   void          AddTrack(Float_t random);
+   Track        *AddTrack(Float_t random);
    void          SetMeasure(UChar_t which, Int_t what);
    void          SetMatrix(UChar_t x, UChar_t y, Float_t what) { if (x<3&&y<3) fMatrix[x][y]=what;}
    void          SetRandomVertex();
@@ -85,6 +138,9 @@ public:
    Float_t       GetTemperature() const { return fTemperature; }
    EventHeader  *GetHeader() { return &fEvtHdr; }
    TClonesArray *GetTracks() const { return fTracks; }
+   TRefArray    *GetHighPt() const {return fHighPt;}
+   TRefArray    *GetMuons()  const {return fMuons;}
+   Track        *GetLastTrack() const {return (Track*)fLastTrack.GetObject();}
    TH1F         *GetHistogram() const { return fH; }
    Int_t         GetMeasure(UChar_t which) { return (which<10)?fMeasures[which]:0; }
    Float_t       GetMatrix(UChar_t x, UChar_t y) { return (x<4&&y<4)?fMatrix[x][y]:0; }
@@ -93,56 +149,6 @@ public:
 };
 
 
-class Track : public TObject {
-
-private:
-   Float_t      fPx;           //X component of the momentum
-   Float_t      fPy;           //Y component of the momentum
-   Float_t      fPz;           //Z component of the momentum
-   Float_t      fRandom;       //A random track quantity
-   Float_t      fMass2;        //The mass square of this particle
-   Float_t      fBx;           //X intercept at the vertex
-   Float_t      fBy;           //Y intercept at the vertex
-   Float_t      fMeanCharge;   //Mean charge deposition of all hits of this track
-   Float_t      fXfirst;       //X coordinate of the first point
-   Float_t      fXlast;        //X coordinate of the last point
-   Float_t      fYfirst;       //Y coordinate of the first point
-   Float_t      fYlast;        //Y coordinate of the last point
-   Float_t      fZfirst;       //Z coordinate of the first point
-   Float_t      fZlast;        //Z coordinate of the last point
-   Float_t      fCharge;       //Charge of this track
-   Float_t      fVertex[3];    //Track vertex position
-   Int_t        fNpoint;       //Number of points for this track
-   Short_t      fValid;        //Validity criterion
-
-public:
-   Track() { }
-   Track(Float_t random);
-   virtual ~Track() { }
-   Float_t       GetPx() const { return fPx; }
-   Float_t       GetPy() const { return fPy; }
-   Float_t       GetPz() const { return fPz; }
-   Float_t       GetPt() const { return TMath::Sqrt(fPx*fPx + fPy*fPy); }
-   Float_t       GetRandom() const { return fRandom; }
-   Float_t       GetBx() const { return fBx; }
-   Float_t       GetBy() const { return fBy; }
-   Float_t       GetMass2() const { return fMass2; }
-   Float_t       GetMeanCharge() const { return fMeanCharge; }
-   Float_t       GetXfirst() const { return fXfirst; }
-   Float_t       GetXlast()  const { return fXlast; }
-   Float_t       GetYfirst() const { return fYfirst; }
-   Float_t       GetYlast()  const { return fYlast; }
-   Float_t       GetZfirst() const { return fZfirst; }
-   Float_t       GetZlast()  const { return fZlast; }
-   Float_t       GetCharge() const { return fCharge; }
-   Float_t       GetVertex(Int_t i=0) {return fVertex[i];}
-   Int_t         GetNpoint() const { return fNpoint; }
-   Short_t       GetValid()  const { return fValid; }
-   virtual void  SetValid(Int_t valid=1) { fValid = valid; }
-
-   ClassDef(Track,1)  //A track segment
-};
-
 
 class HistogramManager {
 
diff --git a/test/MainEvent.cxx b/test/MainEvent.cxx
index 9b933338881..fcaae84d96c 100644
--- a/test/MainEvent.cxx
+++ b/test/MainEvent.cxx
@@ -1,4 +1,4 @@
-// @(#)root/test:$Name:  $:$Id: MainEvent.cxx,v 1.16 2001/04/24 14:32:46 brun Exp $
+// @(#)root/test:$Name:  $:$Id: MainEvent.cxx,v 1.17 2001/07/09 20:38:07 brun Exp $
 // Author: Rene Brun   19/01/97
 
 ////////////////////////////////////////////////////////////////////////
@@ -19,10 +19,13 @@
 //  The statement ***tree->Branch("event", event, 64000,split);*** below
 //  will parse the structure described in Event.h and will make
 //  a new branch for each data member of the class if split is set to 1.
-//    - 5 branches corresponding to the basic types fNtrack,fNseg,fNvertex
-//           ,fFlag and fTemperature.
+//    - 9 branches corresponding to the basic types fType, fNtrack,fNseg,
+//           fNvertex,fFlag,fTemperature,fMeasures,fMatrix,fClosesDistance.
 //    - 3 branches corresponding to the members of the subobject EventHeader.
 //    - one branch for each data member of the class Track of TClonesArray.
+//    - one branch for the TRefArray of high Pt tracks
+//    - one branch for the TRefArray of muon tracks
+//    - one branch for the reference pointer to the last track
 //    - one branch for the object fH (histogram of class TH1F).
 //
 //  if split = 0 only one single branch is created and the complete event
@@ -46,6 +49,9 @@
 //  For each event the event histogram is saved as well as the list
 //  of all tracks.
 //
+//  The two TRefArray contain only references to the original tracks owned by 
+//  the TClonesArray fTracks.
+//
 //  The number of events can be given as the first argument to the program.
 //  By default 400 events are generated.
 //  The compression option can be activated/deactivated via the second argument.
@@ -137,7 +143,9 @@ int main(int argc, char **argv)
    if (arg5 < 100) printev = 1000;
    if (arg5 < 10)  printev = 10000;
 
-   Track::Class()->IgnoreTObjectStreamer();
+   //In this new version of mainEvent, one cannot activate the next statement
+   //because tracks are referenced
+   //Track::Class()->IgnoreTObjectStreamer();
 
 //         Read case
    if (read) {
@@ -236,7 +244,14 @@ int main(int argc, char **argv)
          }
 
          //  Create and Fill the Track objects
-         for (Int_t t = 0; t < ntrack; t++) event->AddTrack(random);
+         Track *track;
+         for (Int_t t = 0; t < ntrack; t++) {
+            track = event->AddTrack(random);
+         
+            //add this track to additional arrays if required
+            if (track->GetPt() > 1) event->GetHighPt()->Add(track);
+            if (track->GetMass2() < 0.11) event->GetMuons()->Add(track);
+         } 
 
          if (write) nb += tree->Fill();  //fill the tree
 
-- 
GitLab