]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCtrackingParamDB.C
DA: Put back the DA rpm description removed by the Doxygen documentation
[u/mrichter/AliRoot.git] / TPC / AliTPCtrackingParamDB.C
1 /// \file AliTPCtrackingParamDB.C
2
3 /****************************************************************************
4  * This macro is used to create a DataBase for the TPC tracking             *
5  * parameterization.                                                        * 
6  * 1) the function CreateAllGeantTracks gives all tracks at the 1st hit of  *
7  *    the TPC                                                               *
8  * 2) the function TrackCompare compares them with track found by the       *
9  *    Kalman filter for the same event and computes efficiency and          *
10  *    resolution on the track parameters for the Kalman filter.             *
11  * 3) the function BuildDataBase calls many functions of AliTPCtrackerParam:*
12  *    - merge results from TrackCompare for many events and compute         *
13  *      average efficiency.                                                 *
14  *    - analyze the pulls of the covariance matrix                          *
15  *    - analyze the dE/dx                                                   *
16  *    - regularize the covariance matrix as a function of the momentum      *
17  *    - write all the informations and the trees with regularized cov.      *
18  *      matrices in the DataBase file.                                      *
19  *                                                                          *
20  *  Origin: A.Dainese, Padova, andrea.dainese@pd,infn.it                    * 
21  *                                                                          *
22  ****************************************************************************/
23
24 #ifndef __CINT__
25 #include "Riostream.h"
26 #include <TFile.h>
27 #include <TTree.h>
28 #include <TStopwatch.h>
29 #include <TObject.h>
30 #include "alles.h"
31 #include "AliRun.h"
32 #include "AliHeader.h"
33 #include "AliGenEventHeader.h"
34 #include "AliMagF.h"
35 #include "AliModule.h"
36 #include "AliArrayI.h"
37 #include "AliDigits.h"
38 #include "AliITS.h"
39 #include "AliTPC.h"
40 #include "AliITSgeom.h"
41 #include "AliITSRecPoint.h"
42 #include "AliITSclusterV2.h"
43 #include "AliITSsimulationFastPoints.h"
44 #include "AliITStrackerV2.h"
45 #include "AliKalmanTrack.h"
46 #include "AliTPCtrackerParam.h"
47 #include "AliTracker.h"
48 #include "AliESD.h"
49 #include "AliRun.h"
50 #include "AliRunLoader.h"
51 #endif
52
53 Int_t TPCParamTracks(const Int_t coll=1,Int_t firstEvent=0,Int_t lastEvent=0);
54 void CreateAllGeantTracks(const Int_t coll=1,Int_t nev=1);
55 void TrackCompare(const Int_t coll,const Double_t Bfield,Int_t n);
56 void BuildDataBase(const Int_t coll,const Double_t Bfield);
57
58 //_____________________________________________________________________________
59 Int_t TPCParamTracks(const Int_t coll,Int_t firstEvent,Int_t lastEvent) {
60 /// Ordinary TPC tracking parameterization
61
62    /**** Initialization of the NewIO *******/
63
64    if (gAlice) {
65       delete AliRunLoader::Instance();
66       delete gAlice; 
67       gAlice=0;
68    }
69
70    AliRunLoader *rl = AliRunLoader::Open("galice.root");
71    if (rl == 0x0) {
72       cerr<<"Can not open session"<<endl;
73       return;
74    }
75    Int_t retval = rl->LoadgAlice();
76    if (retval) {
77       cerr<<"LoadgAlice returned error"<<endl;
78       delete rl;
79       return;
80    }
81    retval = rl->LoadHeader();
82    if (retval) {
83       cerr<<"LoadHeader returned error"<<endl;
84       delete rl;
85       return;
86    }
87    gAlice=rl->GetAliRun();
88        
89
90    TDatabasePDG *DataBase = TDatabasePDG::Instance();
91
92    // Get field from galice.root
93    AliMagF *fiel = TGeoGlobalMagField::Instance()->GetField();
94    Double_t fieval=TMath::Abs((Double_t)fiel->SolenoidField()/10.);
95
96    /**** The TPC corner ********************/
97
98    AliTPCtrackerParam tpcTrackerPar(coll,fieval);
99    tpcTrackerPar.Init();
100
101    /***** The TREE is born *****/
102    
103    TTree *esdTree=new TTree("esdTree","Tree with ESD objects");
104    AliESD *event=0;
105    esdTree->Branch("ESD","AliESD",&event);
106    
107    if(firstEvent>rl->GetNumberOfEvents()) firstEvent=rl->GetNumberOfEvents()-1;
108    if(lastEvent>rl->GetNumberOfEvents())  lastEvent=rl->GetNumberOfEvents()-1;
109    cout<<" Number of events: "<<1+lastEvent-firstEvent<<endl;
110    
111    //<----------------------------------The Loop over events begins
112    TStopwatch timer;
113    Int_t trc;
114    for(Int_t i=firstEvent; i<=lastEvent; i++) { 
115      
116      cerr<<" Processing event number : "<<i<<endl;
117      AliESD *event = new AliESD(); 
118      event->SetRunNumber(gAlice->GetRunNumber());
119      event->SetEventNumber(i);
120      event->SetMagneticField(gAlice->Field()->SolenoidField());
121      rl->GetEvent(i);
122
123      if ( (trc=tpcTrackerPar.BuildTPCtracks(event)) ) {
124        printf("exiting tracker with code %d in event %d\n",trc,i);
125        esdTree->Fill(); delete event;
126        continue;
127      }
128
129      esdTree->Fill();
130      delete event;
131
132    }//<-----------------------------------The Loop over events ends here
133    timer.Stop(); timer.Print();
134
135    //        The AliESDs.root is born
136    TFile *ef = TFile::Open("AliESDs.root","RECREATE"); 
137    if (!ef || !ef->IsOpen()) {cerr<<"Can't open AliESDs.root !\n"; return;}
138
139    esdTree->Write(); //Write the TREE and close everything
140    delete esdTree;
141    ef->Close();
142
143    delete rl;
144
145    return;
146 }
147 //_____________________________________________________________________________
148 void CreateAllGeantTracks(const Int_t coll,Int_t nev) {
149 /// Get all tracks at TPC 1st hit w/o selection and smearing
150
151   cerr<<"\n*******************************************************************\n";
152
153   const Char_t *name="CreateAllGeantTracks";
154   cerr<<'\n'<<name<<"...\n";
155   gBenchmark->Start(name);
156
157   TFile *outfile=TFile::Open(outname,"recreate");
158   TFile *infile =TFile::Open(galice);
159
160   AliTPCtrackerParam tracker(coll,Bfield,n);
161   tracker.AllGeantTracks(); // this is to switch-off selection and smearing
162   tracker.BuildTPCtracks(infile,outfile);
163
164   delete gAlice; gAlice=0;
165
166   infile->Close();
167   outfile->Close();
168
169   gBenchmark->Stop(name);
170   gBenchmark->Show(name);
171
172   return;
173 }
174 //_____________________________________________________________________________
175 void TrackCompare(const Int_t coll,const Double_t Bfield,Int_t n) {
176 /// Compare Kalman tracks with tracks at TPC 1st hit
177
178   cerr<<"\n*******************************************************************\n";
179
180   const Char_t *name="TrackCompare";
181   cerr<<'\n'<<name<<"...\n";
182   gBenchmark->Start(name);
183
184   AliTPCtrackerParam tracker(coll,Bfield,n);
185   tracker.CompareTPCtracks();
186
187   gBenchmark->Stop(name);
188   gBenchmark->Show(name);
189
190   return;
191 }
192 //_____________________________________________________________________________
193 void BuildDataBase(const Int_t coll,const Double_t Bfield) {
194 ///
195
196   cerr<<"\n*******************************************************************\n";
197
198   AliTPCtrackerParam tracker(coll,Bfield);
199
200   // Merge files with cov. matrix and compute average efficiencies
201   cerr<<"\n   --- Merging Events ---\n\n";
202   tracker.MergeEvents(1,5);
203  
204   // Compute the pulls for pions, kaons, electrons
205   cerr<<"\n   --- Analyzing Pulls ---\n\n";
206   tracker.AnalyzePulls("pulls.root");
207
208   // Draw pulls and efficiencies  
209   tracker.DrawPulls("CovMatrixDB_PbPb6000_B0.4T.root",211,0);
210   tracker.DrawEffs("CovMatrixDB_PbPb6000_B0.4T.root",13);
211
212   // Regularize the covariance matrix
213   tracker.RegularizeCovMatrix("regPi.root",211);
214
215   // Analyze the dE/dx
216   tracker.AnalyzedEdx("dEdxPi.root",211);
217
218
219   // Put everything together and create the DB file
220   tracker.MakeDataBase();
221
222   return;
223 }
224
225
226
227
228
229
230