]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtest.C
First implementation of ESD classes (Yu.Belikov)
[u/mrichter/AliRoot.git] / STEER / AliESDtest.C
1 //********************************************************************
2 //     Example of the reconstruction that generates the ESD
3 // Input files: 
4 //   a) AliTPCclusters.root containing the TPC clusters
5 //      (the AliTPCFindClusters.C macro can be used to generate it)
6 //   b) AliITSclustersV2.root containing the ITS clusters
7 //      (the AliITSFindClustersV2.C macro can be used to generate it)
8 // Ouput file:
9 //      AliESDs.root containing the ESD events 
10 //
11 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
12 //********************************************************************
13
14 #ifndef __CINT__
15   #include <Riostream.h>
16   #include "TFile.h"
17   #include "TStopwatch.h"
18
19   #include "AliESD.h"
20   #include "AliTPCParam.h"
21   #include "AliTPCtracker.h"
22   #include "AliITSgeom.h"
23   #include "AliITStrackerV2.h"
24 #endif
25
26 Int_t AliESDtest(Int_t nev=1) { 
27    //File with the TPC clusters
28    TFile *tpccf=TFile::Open("AliTPCclusters.root");
29    if (!tpccf->IsOpen()) {
30       cerr<<"Can't open AliTPCclusters.root !\n"; 
31       return 2;
32    }
33    AliTPCParam *par=(AliTPCParam*)tpccf->Get("75x40_100x60_150x60");
34    if (!par) {cerr<<"Can't get TPC parameters !\n"; return 3;}
35
36    //An instance of the TPC tracker
37    AliTPCtracker tpcTracker(par);
38
39
40    //File with the ITS clusters
41    TFile *itscf=TFile::Open("AliITSclustersV2.root");
42    if (!itscf->IsOpen()) {
43       cerr<<"Can't open AliITSclustersV2.root !\n"; 
44       return 4;
45    }
46    AliITSgeom *geom=(AliITSgeom*)itscf->Get("AliITSgeom");
47    if (!geom) {cerr<<"Can't get AliITSgeom !\n"; return 5;}
48
49    //An instance of the ITS tracker
50    AliITStrackerV2 itsTracker(geom);
51
52    TFile *ef=TFile::Open("AliESDs.root","new");
53    if (!ef->IsOpen()) {cerr<<"Can't AliESDs.root !\n"; return 1;}
54
55    TStopwatch timer;
56    Int_t rc=0;
57    //The loop over events
58    for (Int_t i=0; i<nev; i++) {
59      cerr<<"\n\nProcessing event number : "<<i<<endl;
60      AliESD *event=new AliESD(); 
61  
62      tpcTracker.SetEventNumber(i); 
63      tpcTracker.LoadClusters(tpccf);
64
65      itsTracker.SetEventNumber(i); 
66      itsTracker.LoadClusters(itscf);
67
68      rc+=tpcTracker.Clusters2Tracks(event);
69
70      rc+=itsTracker.Clusters2Tracks(event);
71
72      rc+=itsTracker.PropagateBack(event); 
73      itsTracker.UnloadClusters();
74      
75      rc+=tpcTracker.PropagateBack(event);
76      tpcTracker.UnloadClusters();
77
78      if (rc==0) {
79         Char_t ename[100]; 
80         sprintf(ename,"%d",i);
81         if (!event->Write(ename)) rc++;
82      } 
83      if (rc) {
84         cerr<<"Something bad happened...\n";
85      }
86      delete event;
87    }
88    timer.Stop(); timer.Print();
89
90    delete geom;
91    itscf->Close();
92    delete par;
93    tpccf->Close();
94    ef->Close();
95
96    return rc;
97 }