1 /****************************************************************************
2 * This macro is supposed to do reconstruction in the barrel ALICE trackers *
3 * (Make sure you have TPC digits and ITS hits before using this macro !!!) *
4 * It does the following steps (April 12, 2001): *
5 * 1) TPC cluster finding *
6 * 2) TPC track finding *
7 * 3) ITS cluster finding V2 (via fast points !) *
8 * 4) ITS track finding V2 *
9 * 5) ITS back track propagation V2 *
10 * 6) TPC back track propagation *
11 * (Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch) *
12 ****************************************************************************/
17 #include "AliTPCtracker.h"
20 #include "AliITSgeom.h"
21 #include "AliITSRecPoint.h"
22 #include "AliITSclusterV2.h"
23 #include "AliITSsimulationFastPoints.h"
24 #include "AliITStrackerV2.h"
27 Int_t TPCFindClusters(const Char_t *inname, const Char_t *outname);
28 Int_t TPCFindTracks(const Char_t *inname, const Char_t *outname);
29 Int_t TPCSortTracks(const Char_t *inname, const Char_t *outname);
30 Int_t TPCPropagateBack(const Char_t *inname, const Char_t *outname);
32 Int_t ITSFindClusters(const Char_t *inname, const Char_t *outname);
33 Int_t ITSFindTracks(const Char_t *inname, const Char_t *outname);
34 Int_t ITSPropagateBack(const Char_t *inname, const Char_t *outname);
37 Int_t AliTPCtracknew() {
38 const Char_t *TPCdigName="galice.root";
39 const Char_t *TPCclsName="AliTPCclusters.root";
40 const Char_t *TPCtrkName="AliTPCtracks.root";
42 const Char_t *ITSdigName="galice.root";
43 const Char_t *ITSclsName="AliITSclustersV2.root";
44 const Char_t *ITStrkName="AliITStracksV2.root";
46 AliKalmanTrack::SetConvConst(100/0.299792458/0.2/gAlice->Field()->Factor());
47 // AliKalmanTrack::SetConvConst(100/0.299792458/0.2/2);
48 // ********** Find TPC clusters *********** //
49 if (TPCFindClusters(TPCdigName,TPCclsName)) {
50 cerr<<"Failed to get TPC clusters !\n";
54 // ********** Find TPC tracks *********** //
55 if (TPCFindTracks(TPCclsName,TPCtrkName)) {
56 cerr<<"Failed to get TPC tracks !\n";
60 // ********** Sort and label TPC tracks *********** //
61 if (TPCSortTracks(TPCclsName,TPCtrkName)) {
62 cerr<<"Failed to sort TPC tracks !\n";
66 // ********** Find ITS clusters *********** //
67 if (ITSFindClusters(ITSdigName,ITSclsName)) {
68 cerr<<"Failed to get ITS clusters !\n";
72 // ********** Find ITS tracks *********** //
73 {TFile *clsFile=TFile::Open(ITSclsName);
74 if (ITSFindTracks(TPCtrkName,ITStrkName)) {
75 cerr<<"Failed to get ITS tracks !\n";
80 // ********** Back propagation of the ITS tracks *********** //
81 {TFile *clsFile=TFile::Open(ITSclsName);
82 if (ITSPropagateBack(ITStrkName,TPCtrkName)) {
83 cerr<<"Failed to propagate back the ITS tracks !\n";
89 // ********** Back propagation of the TPC tracks *********** //
90 {TFile *clsFile=TFile::Open(TPCclsName);
91 if (TPCPropagateBack(TPCtrkName,TPCtrkName)) {
92 cerr<<"Failed to propagate back the TPC tracks !\n";
101 Int_t TPCFindClusters(const Char_t *inname, const Char_t *outname) {
103 const Char_t *name="TPCFindClusters";
104 cerr<<'\n'<<name<<"...\n";
105 gBenchmark->Start(name);
107 TFile *out=TFile::Open(outname,"recreate");
108 TFile *in =TFile::Open(inname);
110 AliTPCParam *param=(AliTPCParam *)in->Get("75x40_100x60");
111 if (!param) {cerr<<"TPC parameters have not been found !\n"; return 1;}
115 tpc.Digits2Clusters(out);
119 gBenchmark->Stop(name);
120 gBenchmark->Show(name);
125 Int_t TPCFindTracks(const Char_t *inname, const Char_t *outname) {
127 const Char_t *name="TPCFindTracks";
128 cerr<<'\n'<<name<<"...\n";
129 gBenchmark->Start(name);
130 TFile *out=TFile::Open(outname,"recreate");
131 TFile *in =TFile::Open(inname);
132 AliTPCParam *param=(AliTPCParam *)in->Get("75x40_100x60");
133 if (!param) {cerr<<"TPC parameters have not been found !\n"; return 1;}
135 AliTPCtracker *tracker=new AliTPCtracker(param);
136 rc=tracker->Clusters2Tracks(0,out);
141 gBenchmark->Stop(name);
142 gBenchmark->Show(name);
147 Int_t TPCSortTracks(const Char_t *inname, const Char_t *outname) {
149 const Char_t *name="TPCSortTracks";
150 cerr<<'\n'<<name<<"...\n";
151 gBenchmark->Start(name);
153 TFile *out=TFile::Open(outname);
154 TFile *in =TFile::Open(inname);
155 AliTPCParam *param=(AliTPCParam *)in->Get("75x40_100x60");
156 if (!param) {cerr<<"TPC parameters have not been found !\n"; return 1;}
158 TObjArray tarray(10000);
159 AliTPCtrack *iotrack=0;
163 TTree *tracktree=(TTree*)out->Get("TPCf");
164 TBranch *tbranch=tracktree->GetBranch("tracks");
165 Int_t nentr=(Int_t)tracktree->GetEntries();
166 for (i=0; i<nentr; i++) {
167 iotrack=new AliTPCtrack;
168 tbranch->SetAddress(&iotrack);
169 tracktree->GetEvent(i);
170 tarray.AddLast(iotrack);
175 //assign thacks GEANT labels
177 AliTPCtracker *tracker = new AliTPCtracker(param);
178 tracker->LoadInnerSectors();
179 tracker->LoadOuterSectors();
180 for (i=0; i<nentr; i++) {
181 iotrack=(AliTPCtrack*)tarray.UncheckedAt(i);
182 tracker->CookLabel(iotrack,0.1);
186 //end of GEANT label assignment
188 out=TFile::Open(outname,"recreate");
189 tracktree=new TTree("TPCf","Tree with TPC tracks");
190 tracktree->Branch("tracks","AliTPCtrack",&iotrack,32000,0);
191 for (i=0; i<nentr; i++) {
192 iotrack=(AliTPCtrack*)tarray.UncheckedAt(i);
198 gBenchmark->Stop(name);
199 gBenchmark->Show(name);
204 Int_t ITSFindClusters(const Char_t *inname, const Char_t *outname) {
206 const Char_t *name="ITSFindClusters";
207 cerr<<'\n'<<name<<"...\n";
208 gBenchmark->Start(name);
209 TFile *out=TFile::Open(outname,"recreate");
210 TFile *in =TFile::Open(inname,"update");
212 if (!(gAlice=(AliRun*)in->Get("gAlice"))) {
213 cerr<<"Can't get gAlice !\n";
217 gAlice->GetEvent(ev);
218 AliITS *ITS = (AliITS*)gAlice->GetModule("ITS");
219 if (!ITS) { cerr<<"Can't get the ITS !\n"; return 1;}
221 gAlice->MakeTree("R"); ITS->MakeBranch("R",0);
222 //////////////// Taken from ITSHitsToFastPoints.C ///////////////////////
223 AliITSsimulationFastPoints *sim = new AliITSsimulationFastPoints();
224 for (Int_t i=0;i<3;i++) { ITS->SetSimulationModel(i,sim); }
227 Int_t bgr_ev=Int_t(ev/nsignal);
228 ITS->HitsToFastRecPoints(ev,bgr_ev,size," ","All"," ");
229 //////////////////////////////////////////////////////////////////////////
231 gAlice->GetEvent(ev);
233 AliITSgeom *geom=ITS->GetITSgeom();
237 TClonesArray *clusters=new TClonesArray("AliITSclusterV2",10000);
238 TTree *cTree=new TTree("cTree","ITS clusters");
239 cTree->Branch("Clusters",&clusters);
241 TTree *pTree=gAlice->TreeR();
242 if (!pTree) { cerr<<"Can't get TreeR !\n"; return 1; }
243 TBranch *branch=pTree->GetBranch("ITSRecPoints");
244 if (!branch) { cerr<<"Can't get ITSRecPoints branch !\n"; return 1;}
245 TClonesArray *points=new TClonesArray("AliITSRecPoint",10000);
246 branch->SetAddress(&points);
248 TClonesArray &cl=*clusters;
250 Int_t nentr=(Int_t)pTree->GetEntries();
251 for (Int_t i=0; i<nentr; i++) {
252 if (!pTree->GetEvent(i)) {cTree->Fill(); continue;}
253 Int_t lay,lad,det; geom->GetModuleId(i,lay,lad,det);
254 Float_t x,y,zshift; geom->GetTrans(lay,lad,det,x,y,zshift);
255 Double_t rot[9]; geom->GetRotMatrix(lay,lad,det,rot);
256 Double_t yshift = x*rot[0] + y*rot[1];
257 Int_t ndet=(lad-1)*geom->GetNdetectors(lay) + (det-1);
258 Int_t ncl=points->GetEntriesFast();
260 for (Int_t j=0; j<ncl; j++) {
261 AliITSRecPoint *p=(AliITSRecPoint*)points->UncheckedAt(j);
263 lp[0]=-p->GetX()-yshift; if (lay==1) lp[0]=-lp[0];
264 lp[1]=p->GetZ()+zshift;
265 lp[2]=p->GetSigmaX2();
266 lp[3]=p->GetSigmaZ2();
269 lab[0]=p->GetLabel(0);lab[1]=p->GetLabel(1);lab[2]=p->GetLabel(2);
273 TParticle *part=(TParticle*)gAlice->Particle(label);
275 while (part->P() < 0.005) {
276 Int_t m=part->GetFirstMother();
277 if (m<0) {cerr<<"Primary momentum: "<<part->P()<<endl; break;}
279 part=(TParticle*)gAlice->Particle(label);
281 if (lab[1]<0) lab[1]=label;
282 else if (lab[2]<0) lab[2]=label;
283 else cerr<<"No empty labels !\n";
285 new(cl[j]) AliITSclusterV2(lab,lp);
287 cTree->Fill(); clusters->Delete();
291 cerr<<"Number of clusters: "<<nclusters<<endl;
293 delete cTree; delete clusters; delete points;
295 delete gAlice; gAlice=0;
298 gBenchmark->Stop(name);
299 gBenchmark->Show(name);
304 Int_t ITSFindTracks(const Char_t *inname, const Char_t *outname) {
306 const Char_t *name="ITSFindTracks";
307 cerr<<'\n'<<name<<"...\n";
308 gBenchmark->Start(name);
310 AliITSgeom *geom=(AliITSgeom*)gFile->Get("AliITSgeom");
311 if (!geom) { cerr<<"can't get ITS geometry !\n"; return 1;}
312 AliITStrackerV2 tracker(geom);
314 TFile *out=TFile::Open(outname,"recreate");
315 TFile *in =TFile::Open(inname);
316 rc=tracker.Clusters2Tracks(in,out);
320 gBenchmark->Stop(name);
321 gBenchmark->Show(name);
326 Int_t ITSPropagateBack(const Char_t *inname, const Char_t *outname) {
328 const Char_t *name="ITSPropagateBack";
329 cerr<<'\n'<<name<<"...\n";
330 gBenchmark->Start(name);
332 AliITSgeom *geom=(AliITSgeom*)gFile->Get("AliITSgeom");
333 if (!geom) { cerr<<"can't get ITS geometry !\n"; return 1;}
334 AliITStrackerV2 tracker(geom);
336 TFile *out=TFile::Open(outname,"update");
337 TFile *in =TFile::Open(inname);
338 rc=tracker.PropagateBack(in,out);
342 gBenchmark->Stop(name);
343 gBenchmark->Show(name);
348 Int_t TPCPropagateBack(const Char_t *inname, const Char_t *outname) {
350 const Char_t *name="TPCPropagateBack";
351 cerr<<'\n'<<name<<"...\n";
352 gBenchmark->Start(name);
354 AliTPCParam *param=(AliTPCParam *)gFile->Get("75x40_100x60");
355 if (!param) {cerr<<"TPC parameters have not been found !\n"; return 1;}
356 AliTPCtracker *tracker=new AliTPCtracker(param);
358 TFile *out=TFile::Open(outname,"update");
359 TFile *in =TFile::Open(inname);
360 rc=tracker->PropagateBack(in,out);
365 gBenchmark->Stop(name);
366 gBenchmark->Show(name);