]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ReadESDfriend.C
Adding the new MUONcalib library (Ivana)
[u/mrichter/AliRoot.git] / STEER / ReadESDfriend.C
1 //********************************************************************
2 //  Example of accessing the information stored in ESD friends
3 //      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
4 //********************************************************************
5
6 #if !defined( __CINT__) || defined(__MAKECINT__)
7   #include <Riostream.h>
8   #include <TFile.h>
9   #include <TChain.h>
10
11   #include "AliESD.h"
12   #include "AliESDfriend.h"
13   #include "AliTrackPointArray.h"
14 #endif
15
16 void ReadESDfriend(Bool_t readFriend=kTRUE) {
17    Char_t *name[]={
18      //Put here the names of the ESD files to be chained
19      "AliESDs.root"
20    };
21    Int_t n=sizeof(name)/sizeof(Char_t *); 
22    TChain *esdTree=new TChain("esdTree");
23    for (Int_t i=0; i<n; i++) esdTree->AddFile(name[i]);
24
25    AliESD *ev=0;
26    esdTree->SetBranchAddress("ESD",&ev);
27
28    // Attach the branch with ESD friends
29    AliESDfriend *evf=0;
30    if (readFriend) {
31       esdTree->SetBranchStatus("ESDfriend*",1);
32       esdTree->SetBranchAddress("ESDfriend.",&evf);
33    }
34
35    Int_t nev=esdTree->GetEntries();
36    for (Int_t i=0; i<nev; i++) {
37        esdTree->GetEntry(i);
38         
39        cout<<endl<<"Event number: "<<i<<endl;
40        Int_t n=ev->GetNumberOfTracks();
41        cout<<"Number of tracks: "<<n<<endl;
42
43        ev->SetESDfriend(evf); //Attach the friend to the ESD
44
45     // Now the attached information can be accessed via pointer to ESD.
46     // Example: indices of the TPC clusters associated with the track number 0.
47        if (n > 0) {
48           const AliESDtrack *t=ev->GetTrack(0);
49           Int_t idx[AliESDfriendTrack::kMaxTPCcluster]; 
50           n=t->GetTPCclusters(idx);
51           cout<<"Track number 0"<<endl;
52           cout<<"   Number of TPC clusters: "<<n<<endl;
53           cout<<"   Index of the 7th TPC cluster: "<<idx[7]<<endl;
54           UChar_t map=t->GetITSClusterMap();
55           cout<<"   ITS cluster map (from SPDs to SSDs): ";
56           for (Int_t i=0; i<6; i++) cout<<TESTBIT(map,i)<<' ';
57           cout<<endl;
58
59     // Example: track points associated with the track number 0.
60           const AliTrackPointArray *pa=t->GetTrackPointArray();
61           if (pa != 0) {
62              n=pa->GetNPoints();
63              const Float_t *x=pa->GetX();
64              cout<<"   Number of track points: "<<n<<endl;
65              if (n>7)
66              cout<<"   X coordinate of the 7th track point: "<<x[7]<<endl;
67           }
68        }
69
70        delete ev;  ev=0;
71        delete evf; evf=0;
72
73    }
74
75    delete esdTree;
76 }