]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSTrackSegmentMakerv1.cxx
Remove AliTRDv2
[u/mrichter/AliRoot.git] / PHOS / AliPHOSTrackSegmentMakerv1.cxx
CommitLineData
d15a28e7 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
b2a60966 15/* $Id$ */
702ab87e 16
17/* History of cvs commits:
18 *
19 * $Log$
15e85efa 20 * Revision 1.78 2005/11/18 13:04:51 hristov
21 * Bug fix
22 *
a494460a 23 * Revision 1.77 2005/11/17 23:34:36 hristov
24 * Corrected logics
25 *
f74a631d 26 * Revision 1.76 2005/11/17 22:29:12 hristov
27 * Faster version, no attempt to match tracks outside the PHOS acceptance
28 *
0d9aa319 29 * Revision 1.75 2005/11/17 12:35:27 hristov
30 * Use references instead of objects. Avoid to create objects when they are not really needed
31 *
7b51037f 32 * Revision 1.74 2005/07/08 14:01:36 hristov
33 * Tracking in non-uniform nmagnetic field (Yu.Belikov)
34 *
c84a5e9e 35 * Revision 1.73 2005/05/28 14:19:05 schutz
36 * Compilation warnings fixed by T.P.
37 *
702ab87e 38 */
39
d15a28e7 40//_________________________________________________________________________
b2a60966 41// Implementation version 1 of algorithm class to construct PHOS track segments
f035f6ce 42// Track segment for PHOS is list of
194f9939 43// EMC RecPoint + (possibly) CPV RecPoint
a4e98857 44// To find TrackSegments we do the following:
45// for each EMC RecPoints we look at
194f9939 46// CPV RecPoints in the radious fRcpv.
a4e98857 47// If there is such a CPV RecPoint,
194f9939 48// we make "Link" it is just indexes of EMC and CPV RecPoint and distance
a4e98857 49// between them in the PHOS plane.
50// Then we sort "Links" and starting from the
51// least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
52// new TrackSegment.
194f9939 53// If there is no CPV RecPoint we make TrackSegment
a4e98857 54// consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
f444a19f 55//// In principle this class should be called from AliPHOSReconstructor, but
a4e98857 56// one can use it as well in standalone mode.
57// Use case:
fc12304f 58// root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
a4e98857 59// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
fc12304f 60// // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
61// // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
a4e98857 62// root [1] t->ExecuteTask()
a4e98857 63// root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
64// root [4] t->ExecuteTask("deb all time")
f035f6ce 65//
fc12304f 66//*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
b2a60966 67//
d15a28e7 68
69// --- ROOT system ---
2731cd1e 70#include "TTree.h"
2731cd1e 71#include "TBenchmark.h"
d15a28e7 72
21cd0c07 73// --- Standard library ---
194f9939 74#include "Riostream.h"
d15a28e7 75// --- AliRoot header files ---
e957fea8 76#include "AliPHOSGeometry.h"
d15a28e7 77#include "AliPHOSTrackSegmentMakerv1.h"
78#include "AliPHOSTrackSegment.h"
79#include "AliPHOSLink.h"
7b7c1533 80#include "AliPHOSGetter.h"
aa0b9641 81#include "AliESD.h"
82#include "AliESDtrack.h"
d15a28e7 83
84ClassImp( AliPHOSTrackSegmentMakerv1)
85
86
87//____________________________________________________________________________
2bd5457f 88 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
d15a28e7 89{
7b7c1533 90 // default ctor (to be used mainly by Streamer)
91
8d0f3f77 92 InitParameters() ;
92f521a9 93 fDefaultInit = kTRUE ;
d15a28e7 94}
7b7c1533 95
9f616d61 96//____________________________________________________________________________
88cb7938 97 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString alirunFileName, const TString eventFolderName)
98 :AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName)
2731cd1e 99{
100 // ctor
2731cd1e 101
8d0f3f77 102 InitParameters() ;
7b7c1533 103 Init() ;
92f521a9 104 fDefaultInit = kFALSE ;
aa0b9641 105 fESD = 0;
2731cd1e 106}
98cbd830 107
2731cd1e 108//____________________________________________________________________________
109 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
110{
111 // dtor
92f521a9 112 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
88cb7938 113 if (!fDefaultInit)
114 delete fLinkUpArray ;
d15a28e7 115}
9f616d61 116
8d0f3f77 117
fc12304f 118//____________________________________________________________________________
119const TString AliPHOSTrackSegmentMakerv1::BranchName() const
120{
88cb7938 121
122 return GetName() ;
fc12304f 123}
124
d15a28e7 125//____________________________________________________________________________
2731cd1e 126void AliPHOSTrackSegmentMakerv1::FillOneModule()
9f616d61 127{
f035f6ce 128 // Finds first and last indexes between which
129 // clusters from one PHOS module are
88cb7938 130
131 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
fc12304f 132
fbf811ec 133 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
134 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
9688c1dd 135
2731cd1e 136 //First EMC clusters
7b7c1533 137 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
2731cd1e 138 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
29b077b5 139 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
2731cd1e 140 fEmcLast ++) ;
141
2731cd1e 142 //Now CPV clusters
7b7c1533 143 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
6ad0bfa0 144
2731cd1e 145 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
88cb7938 146 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
147 fCpvLast ++) ;
9688c1dd 148
d15a28e7 149}
7b7c1533 150
d15a28e7 151//____________________________________________________________________________
194f9939 152Float_t AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSCpvRecPoint * cpvClu, Int_t &trackindex) const
d15a28e7 153{
194f9939 154 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
a4e98857 155 // Clusters are sorted in "rows" and "columns" of width 1 cm
f035f6ce 156
743cb288 157 //Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
2731cd1e 158 // if you change this value, change it as well in xxxRecPoint::Compare()
194f9939 159 Float_t distance2Cpv = fRcpv ;
743cb288 160 Float_t distance2Track = fRtpc ;
194f9939 161
162 trackindex = -1 ; // closest track within fRCpv
163
aa0b9641 164 TVector3 vecEmc ; // Local position of EMC recpoint
165 TVector3 vecCpv ; // Local position of CPV recpoint propagated to EMC
166 TVector3 vecDist ; // Distance between local positions of two points
2731cd1e 167
168 emcClu->GetLocalPosition(vecEmc) ;
bfc17d18 169 cpvClu->GetLocalPosition(vecCpv) ;
2731cd1e 170
194f9939 171 //toofar = kTRUE ;
2731cd1e 172 if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){
2731cd1e 173
194f9939 174 // Find EMC-CPV distance
175 distance2Cpv = (vecCpv - vecEmc).Mag() ;
176
aa0b9641 177 if (fESD != 0x0) {
0d9aa319 178 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
179 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
180
181 Double_t rPHOS = geom->GetIPtoCrystalSurface();
182
183 //PH Acceptance boundaries for each PHOS module
184 Int_t nModules = geom->GetNModules();
185 Double_t * thmin = new Double_t[nModules];// theta min
186 Double_t * thmax = new Double_t[nModules];// theta max
187 Double_t * phmin = new Double_t[nModules];// phi min
188 Double_t * phmax = new Double_t[nModules];// phi max
189
190 for (Int_t imod=0; imod<nModules; imod++) {
191 geom->EmcModuleCoverage(imod,
192 thmin[imod],thmax[imod],
193 phmin[imod],phmax[imod]);
194 }
195
aa0b9641 196 // Extrapolate the global track direction if any to CPV and find the closest track
197 Int_t nTracks = fESD->GetNumberOfTracks();
198 Int_t iClosestTrack = -1;
aa0b9641 199 Double_t minDistance = 1e6;
200 Double_t pxyz[3], xyz[3];
0d9aa319 201 TVector3 inPHOS; //PH Used to calculate theta and phi
202
203 //PH Loop on tracks
aa0b9641 204 AliESDtrack *track;
205 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
206 track = fESD->GetTrack(iTrack);
c84a5e9e 207 if (!track->GetXYZAt(rPHOS, fESD->GetMagneticField(), xyz))
23904d16 208 continue; //track coord on the cylinder of PHOS radius
75408399 209 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
23904d16 210 continue;
0d9aa319 211 //PH Here one has to cut out the tracks which are not inside the PHOS
212 //PH acceptance
213 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
214 Double_t inPhi = inPHOS.Phi();
215 Double_t inTheta = inPHOS.Theta();
216
f74a631d 217 Bool_t skip = kTRUE;
a494460a 218 for (Int_t imod=0; imod<nModules; imod++) {
0d9aa319 219 //PH Loop on modules to check if the track enters in the acceptance
f74a631d 220 if (thmin[imod] < inTheta && thmax[imod] > inTheta &&
221 phmin[imod] < inPhi && phmax[imod] > inPhi) {
222 skip = kFALSE;
0d9aa319 223 break;
224 }
225 }
226 if (skip) continue; //PH Skip, if not in the PHOS acceptance
227
c84a5e9e 228 if (!track->GetPxPyPzAt(rPHOS, fESD->GetMagneticField(), pxyz))
23904d16 229 continue; // track momentum ibid.
7b51037f 230 PropagateToPlane(vecDist,xyz,pxyz,"CPV",cpvClu->GetPHOSMod());
e718c158 231 // Info("GetDistanceInPHOSPlane","Track %d propagation to CPV = (%f,%f,%f)",
232 // iTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
aa0b9641 233 vecDist -= vecCpv;
194f9939 234 distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
aa0b9641 235 // Find the closest track to the EMC recpoint
194f9939 236 if (distance2Track < minDistance) {
237 minDistance = distance2Track;
aa0b9641 238 iClosestTrack = iTrack;
239 }
240 }
241
0d9aa319 242 delete [] thmin;
243 delete [] thmax;
244 delete [] phmin;
245 delete [] phmax;
246
aa0b9641 247 if (iClosestTrack != -1) {
248 track = fESD->GetTrack(iClosestTrack);
c84a5e9e 249 if (track->GetPxPyPzAt(rPHOS, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
aa0b9641 250 TVector3 vecCpvGlobal; // Global position of the CPV recpoint
aa0b9641 251 geom->GetGlobal((AliRecPoint*)cpvClu,vecCpvGlobal);
252 for (Int_t ixyz=0; ixyz<3; ixyz++)
253 xyz[ixyz] = vecCpvGlobal[ixyz];
7b51037f 254 PropagateToPlane(vecDist,xyz,pxyz,"EMC",cpvClu->GetPHOSMod());
bfc17d18 255// Info("GetDistanceInPHOSPlane","Track %d propagation to EMC = (%f,%f,%f)",
256// iClosestTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
aa0b9641 257 vecDist -= vecEmc;
194f9939 258 distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
23904d16 259 }
aa0b9641 260 }
194f9939 261// } else {
262// // If no ESD exists, than simply find EMC-CPV distance
263// distance = (vecCpv - vecEmc).Mag() ;
264
743cb288 265 //if(distance2Track < fRcpv + 2*delta )
266 if(distance2Track < fRtpc )
194f9939 267 trackindex = iClosestTrack ;
268 // toofar = kFALSE ;
aa0b9641 269 }
194f9939 270 // Info("GetDistanceInPHOSPlane","cpv-emc distance is %f cm",
271 // distance);
bfc17d18 272 }
d15a28e7 273
194f9939 274 return distance2Cpv ;
aa0b9641 275}
276
277//____________________________________________________________________________
7b51037f 278void AliPHOSTrackSegmentMakerv1::PropagateToPlane(TVector3& globalIntersection,
279 Double_t *x,
280 Double_t *p,
281 const char *det,
282 Int_t moduleNumber) const
aa0b9641 283{
284 // Propagate a straight-line track from the origin point x
bfc17d18 285 // along the direction p to the CPV or EMC module moduleNumber
aa0b9641 286 // Returns a local position of such a propagation
287
288 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
289 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
7b51037f 290 TVector3 moduleCenter;
291 geom->GetModuleCenter(moduleCenter,det,moduleNumber);
292 TVector3 vertex; vertex.SetXYZ(x[0],x[1],x[2]);
293 TVector3 direction; direction.SetXYZ(p[0],p[1],p[2]);
aa0b9641 294
bfc17d18 295// Info("PropagateToCPV","Center of the %s module %d is (%f,%f,%f)",
296// det,moduleNumber,moduleCenter[0],moduleCenter[1],moduleCenter[2]);
395f4eea 297
aa0b9641 298 Double_t time = (moduleCenter.Mag2() - vertex.Dot(moduleCenter)) /
299 (direction.Dot(moduleCenter));
7b51037f 300 vertex += direction*time;
301 geom->Global2Local(globalIntersection,vertex,moduleNumber);
d15a28e7 302}
303
7b7c1533 304//____________________________________________________________________________
305void AliPHOSTrackSegmentMakerv1::Init()
306{
307 // Make all memory allocations that are not possible in default constructor
308
adcca1e6 309 AliPHOSGetter* gime = AliPHOSGetter::Instance();
310 if(!gime)
311 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
fbf811ec 312
7b7c1533 313 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
88cb7938 314 if ( !gime->TrackSegmentMaker() ) {
315 gime->PostTrackSegmentMaker(this);
316 }
7b7c1533 317}
318
8d0f3f77 319//____________________________________________________________________________
320void AliPHOSTrackSegmentMakerv1::InitParameters()
321{
e957fea8 322 //Initializes parameters
743cb288 323 fRcpv = 10. ;
324 fRtpc = 4. ;
8d0f3f77 325 fEmcFirst = 0 ;
326 fEmcLast = 0 ;
327 fCpvFirst = 0 ;
328 fCpvLast = 0 ;
329 fLinkUpArray = 0 ;
adcca1e6 330 fWrite = kTRUE ;
88cb7938 331 fTrackSegmentsInRun = 0 ;
eabde521 332 SetEventRange(0,-1) ;
8d0f3f77 333}
334
335
d15a28e7 336//____________________________________________________________________________
baef0810 337void AliPHOSTrackSegmentMakerv1::MakeLinks()const
d15a28e7 338{
194f9939 339 // Finds distances (links) between all EMC and CPV clusters,
fbf811ec 340 // which are not further apart from each other than fRcpv
f035f6ce 341 // and sort them in accordance with this distance
9688c1dd 342
88cb7938 343 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
fbf811ec 344 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
345 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
7b7c1533 346
2731cd1e 347 fLinkUpArray->Clear() ;
2731cd1e 348
2bb500e5 349 AliPHOSCpvRecPoint * cpv ;
92862013 350 AliPHOSEmcRecPoint * emcclu ;
28c3a259 351
d15a28e7 352 Int_t iLinkUp = 0 ;
353
28c3a259 354 Int_t iEmcRP;
2731cd1e 355 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
29b077b5 356 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
2731cd1e 357
194f9939 358 //Bool_t toofar ;
2731cd1e 359 Int_t iCpv = 0 ;
360 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
28c3a259 361
2bb500e5 362 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
194f9939 363 Int_t track = -1 ;
364 Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, track) ;
365 // if(toofar)
366 // continue ;
fbf811ec 367 if(r < fRcpv) {
194f9939 368 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(r, iEmcRP, iCpv, track) ;
28c3a259 369 }
d15a28e7 370 }
28c3a259 371 }
d15a28e7 372
9688c1dd 373 fLinkUpArray->Sort() ; //first links with smallest distances
d15a28e7 374}
28c3a259 375
d15a28e7 376//____________________________________________________________________________
2731cd1e 377void AliPHOSTrackSegmentMakerv1::MakePairs()
6ad0bfa0 378{
f035f6ce 379 // Using the previously made list of "links", we found the smallest link - i.e.
a4e98857 380 // link with the least distance between EMC and CPV and pointing to still
f035f6ce 381 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
382 // remove them from the list of "unassigned".
88cb7938 383
384 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
385
fbf811ec 386 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
387 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
88cb7938 388 TClonesArray * trackSegments = gime->TrackSegments();
9688c1dd 389
01a599c9 390 //Make arrays to mark clusters already chosen
2731cd1e 391 Int_t * emcExist = 0;
392 if(fEmcLast > fEmcFirst)
393 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
394
395 Int_t index;
396 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
397 emcExist[index] = 1 ;
398
399 Bool_t * cpvExist = 0;
400 if(fCpvLast > fCpvFirst)
401 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
402 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
403 cpvExist[index] = kTRUE ;
404
2731cd1e 405
406 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
2731cd1e 407 TIter nextUp(fLinkUpArray) ;
d15a28e7 408
d15a28e7 409 AliPHOSLink * linkUp ;
9688c1dd 410
2bb500e5 411 AliPHOSCpvRecPoint * nullpointer = 0 ;
9688c1dd 412
29b077b5 413 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
9688c1dd 414
194f9939 415 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
d15a28e7 416
194f9939 417 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
418 new ((* trackSegments)[fNTrackSegments])
419 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
420 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
421 linkUp->GetTrack()) ;
422
88cb7938 423 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
424 fNTrackSegments++ ;
88cb7938 425 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
426 //mark CPV recpoint as already used
194f9939 427 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
428 } //if CpvUp still exist
28c3a259 429 }
88cb7938 430 }
28c3a259 431
194f9939 432 //look through emc recPoints left without CPV
2731cd1e 433 if(emcExist){ //if there is emc rec point
434 Int_t iEmcRP ;
435 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
436 if(emcExist[iEmcRP] > 0 ){
88cb7938 437 new ((*trackSegments)[fNTrackSegments])
438 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
439 nullpointer) ;
440 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
441 fNTrackSegments++;
2731cd1e 442 }
d15a28e7 443 }
d15a28e7 444 }
780a31c1 445 delete [] emcExist ;
446 delete [] cpvExist ;
d15a28e7 447}
448
449//____________________________________________________________________________
eabde521 450void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
d15a28e7 451{
eabde521 452 // Steering method to perform track segment construction for events
453 // in the range from fFirstEvent to fLastEvent.
454 // This range is optionally set by SetEventRange().
455 // if fLastEvent=-1 (by default), then process events until the end.
88cb7938 456
2731cd1e 457 if(strstr(option,"tim"))
b3f97575 458 gBenchmark->Start("PHOSTSMaker");
7b7c1533 459
460 if(strstr(option,"print")) {
88cb7938 461 Print() ;
7b7c1533 462 return ;
463 }
88cb7938 464
adcca1e6 465 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
fb43ada4 466
55ea5766 467 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
7b7c1533 468
0588c896 469 if (fLastEvent == -1)
470 fLastEvent = gime->MaxEvent() - 1 ;
471 else
472 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
eabde521 473 Int_t nEvents = fLastEvent - fFirstEvent + 1;
01a599c9 474
fb43ada4 475 Int_t ievent ;
476 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
aa0b9641 477 gime->Event(ievent,"R") ;
0588c896 478 //Make some initializations
55ea5766 479 fNTrackSegments = 0 ;
480 fEmcFirst = 0 ;
481 fEmcLast = 0 ;
482 fCpvFirst = 0 ;
483 fCpvLast = 0 ;
88cb7938 484
485 gime->TrackSegments()->Clear();
55ea5766 486
9688c1dd 487 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
488
88cb7938 489 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
2731cd1e 490 FillOneModule() ;
2731cd1e 491 MakeLinks() ;
2731cd1e 492 MakePairs() ;
2731cd1e 493 }
28c3a259 494
90cceaf6 495 WriteTrackSegments() ;
7b7c1533 496
2731cd1e 497 if(strstr(option,"deb"))
88cb7938 498 PrintTrackSegments(option);
94de8339 499
500 //increment the total number of track segments per run
88cb7938 501 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
2731cd1e 502 }
88cb7938 503
2731cd1e 504 if(strstr(option,"tim")){
505 gBenchmark->Stop("PHOSTSMaker");
21cd0c07 506 Info("Exec", "took %f seconds for making TS %f seconds per event",
88cb7938 507 gBenchmark->GetCpuTime("PHOSTSMaker"),
eabde521 508 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
88cb7938 509 }
adcca1e6 510 if(fWrite) //do not unload in "on flight" mode
511 Unload();
88cb7938 512}
513
514//____________________________________________________________________________
515void AliPHOSTrackSegmentMakerv1::Unload()
516{
e957fea8 517 // Unloads the task from the folder
88cb7938 518 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
519 gime->PhosLoader()->UnloadRecPoints() ;
520 gime->PhosLoader()->UnloadTracks() ;
d15a28e7 521}
7b7c1533 522
d15a28e7 523//____________________________________________________________________________
702ab87e 524void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
a4e98857 525{
baef0810 526 // Print TrackSegmentMaker parameters
527
21cd0c07 528 TString message("") ;
7b7c1533 529 if( strcmp(GetName(), "") != 0 ) {
21cd0c07 530 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
531 message += "Making Track segments\n" ;
532 message += "with parameters:\n" ;
194f9939 533 message += " Maximal EMC - CPV distance (cm) %f\n" ;
21cd0c07 534 message += "============================================\n" ;
535 Info("Print", message.Data(),fRcpv) ;
2731cd1e 536 }
537 else
21cd0c07 538 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
d15a28e7 539}
7b7c1533 540
d15a28e7 541//____________________________________________________________________________
90cceaf6 542void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
a4e98857 543{
f035f6ce 544 // Writes found TrackSegments to TreeR. Creates branches
545 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
546 // In the former branch found TrackSegments are stored, while
547 // in the latter all parameters, with which TS were made.
548 // ROOT does not allow overwriting existing branches, therefore
a4e98857 549 // first we check, if branches with the same title already exist.
f035f6ce 550 // If yes - exits without writing.
88cb7938 551
552 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
55ea5766 553
fbf811ec 554 TClonesArray * trackSegments = gime->TrackSegments() ;
55ea5766 555 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
8d0f3f77 556
adcca1e6 557 if(fWrite){ //We write TreeT
558 TTree * treeT = gime->TreeT();
559
560 //First TS
561 Int_t bufferSize = 32000 ;
562 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
563 tsBranch->Fill() ;
564
565 gime->WriteTracks("OVERWRITE");
566 gime->WriteTrackSegmentMaker("OVERWRITE");
567 }
2731cd1e 568}
98cbd830 569
98cbd830 570
2731cd1e 571//____________________________________________________________________________
a4e98857 572void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
573{
f035f6ce 574 // option deb - prints # of found TrackSegments
575 // option deb all - prints as well indexed of found RecParticles assigned to the TS
9688c1dd 576
88cb7938 577 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
21cd0c07 578
88cb7938 579 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
580 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
581 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
582
2731cd1e 583 if(strstr(option,"all")) { // printing found TS
88cb7938 584 printf("TrackSegment # EMC RP# CPV RP#\n") ;
2731cd1e 585 Int_t index;
7b7c1533 586 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
587 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
88cb7938 588 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;
2731cd1e 589 }
d15a28e7 590 }
2731cd1e 591}