]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSTrackSegmentMakerv1.cxx
Added new classes
[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$ */
d15a28e7 16//_________________________________________________________________________
b2a60966 17// Implementation version 1 of algorithm class to construct PHOS track segments
f035f6ce 18// Track segment for PHOS is list of
19// EMC RecPoint + (possibly) CPV RecPoint + (possibly) PPSD RecPoint
a4e98857 20// To find TrackSegments we do the following:
21// for each EMC RecPoints we look at
fbf811ec 22// CPV/PPSD RecPoints in the radious fRcpv.
a4e98857 23// If there is such a CPV RecPoint,
24// we make "Link" it is just indexes of EMC and CPV/PPSD RecPoint and distance
25// between them in the PHOS plane.
26// Then we sort "Links" and starting from the
27// least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
28// new TrackSegment.
29// If there is no CPV/PPSD RecPoint we make TrackSegment
30// consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
21cd0c07 31//// In principle this class should be called from AliPHOSReconstructioner, but
a4e98857 32// one can use it as well in standalone mode.
33// Use case:
fc12304f 34// root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
a4e98857 35// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
fc12304f 36// // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
37// // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
a4e98857 38// root [1] t->ExecuteTask()
39// root [2] t->SetMaxEmcPpsdDistance(5)
40// root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
41// root [4] t->ExecuteTask("deb all time")
f035f6ce 42//
fc12304f 43//*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
b2a60966 44//
d15a28e7 45
46// --- ROOT system ---
2731cd1e 47#include "TTree.h"
2731cd1e 48#include "TBenchmark.h"
d15a28e7 49
21cd0c07 50// --- Standard library ---
d15a28e7 51
52// --- AliRoot header files ---
e957fea8 53#include "AliPHOSGeometry.h"
d15a28e7 54#include "AliPHOSTrackSegmentMakerv1.h"
55#include "AliPHOSTrackSegment.h"
56#include "AliPHOSLink.h"
7b7c1533 57#include "AliPHOSGetter.h"
aa0b9641 58#include "AliESD.h"
59#include "AliESDtrack.h"
d15a28e7 60
61ClassImp( AliPHOSTrackSegmentMakerv1)
62
63
64//____________________________________________________________________________
2bd5457f 65 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
d15a28e7 66{
7b7c1533 67 // default ctor (to be used mainly by Streamer)
68
8d0f3f77 69 InitParameters() ;
92f521a9 70 fDefaultInit = kTRUE ;
d15a28e7 71}
7b7c1533 72
9f616d61 73//____________________________________________________________________________
88cb7938 74 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString alirunFileName, const TString eventFolderName)
75 :AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName)
2731cd1e 76{
77 // ctor
2731cd1e 78
8d0f3f77 79 InitParameters() ;
7b7c1533 80 Init() ;
92f521a9 81 fDefaultInit = kFALSE ;
aa0b9641 82 fESD = 0;
2731cd1e 83}
98cbd830 84
2731cd1e 85//____________________________________________________________________________
86 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
87{
88 // dtor
92f521a9 89 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
88cb7938 90 if (!fDefaultInit)
91 delete fLinkUpArray ;
d15a28e7 92}
9f616d61 93
8d0f3f77 94
fc12304f 95//____________________________________________________________________________
96const TString AliPHOSTrackSegmentMakerv1::BranchName() const
97{
88cb7938 98
99 return GetName() ;
fc12304f 100}
101
d15a28e7 102//____________________________________________________________________________
2731cd1e 103void AliPHOSTrackSegmentMakerv1::FillOneModule()
9f616d61 104{
f035f6ce 105 // Finds first and last indexes between which
106 // clusters from one PHOS module are
88cb7938 107
108 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
fc12304f 109
fbf811ec 110 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
111 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
9688c1dd 112
2731cd1e 113 //First EMC clusters
7b7c1533 114 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
2731cd1e 115 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
29b077b5 116 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
2731cd1e 117 fEmcLast ++) ;
118
2731cd1e 119 //Now CPV clusters
7b7c1533 120 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
6ad0bfa0 121
2731cd1e 122 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
88cb7938 123 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
124 fCpvLast ++) ;
9688c1dd 125
d15a28e7 126}
7b7c1533 127
d15a28e7 128//____________________________________________________________________________
aa0b9641 129Float_t AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSRecPoint * cpvClu, Bool_t &toofar) const
d15a28e7 130{
b2a60966 131 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
a4e98857 132 // Clusters are sorted in "rows" and "columns" of width 1 cm
f035f6ce 133
2731cd1e 134 Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
135 // if you change this value, change it as well in xxxRecPoint::Compare()
aa0b9641 136 Float_t distance = fRcpv ;
d15a28e7 137
aa0b9641 138 TVector3 vecEmc ; // Local position of EMC recpoint
139 TVector3 vecCpv ; // Local position of CPV recpoint propagated to EMC
140 TVector3 vecDist ; // Distance between local positions of two points
2731cd1e 141
142 emcClu->GetLocalPosition(vecEmc) ;
143 cpvClu->GetLocalPosition(vecCpv) ;
144
aa0b9641 145 toofar = kTRUE ;
2731cd1e 146 if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){
2731cd1e 147
aa0b9641 148 if (fESD != 0x0) {
149 // Extrapolate the global track direction if any to CPV and find the closest track
150 Int_t nTracks = fESD->GetNumberOfTracks();
151 Int_t iClosestTrack = -1;
152 Double_t distance;
153 Double_t minDistance = 1e6;
154 Double_t pxyz[3], xyz[3];
155 AliESDtrack *track;
156 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
157 track = fESD->GetTrack(iTrack);
158 track->GetOuterXYZ(xyz); // track coord on the cylinder of PHOS radius
75408399 159 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
160 continue;
aa0b9641 161 track->GetOuterPxPyPz(pxyz); // track momentum ibid.
162 vecDist = PropagateToCPV(xyz,pxyz,cpvClu->GetPHOSMod());
163 vecDist -= vecCpv;
164 distance = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
165 // Find the closest track to the EMC recpoint
166 if (distance < minDistance) {
167 minDistance = distance;
168 iClosestTrack = iTrack;
169 }
170 }
171
172 if (iClosestTrack != -1) {
173 track = fESD->GetTrack(iClosestTrack);
174 track->GetOuterPxPyPz(pxyz); // track momentum ibid.
175 TVector3 vecCpvGlobal; // Global position of the CPV recpoint
176 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
177 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
178 geom->GetGlobal((AliRecPoint*)cpvClu,vecCpvGlobal);
179 for (Int_t ixyz=0; ixyz<3; ixyz++)
180 xyz[ixyz] = vecCpvGlobal[ixyz];
181 vecDist = PropagateToCPV(xyz,pxyz,cpvClu->GetPHOSMod());
182 vecDist -= vecEmc;
183 distance = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
184 }
185 } else {
186 // If no ESD, than simply find EMC-CPV distance
187 distance = (vecCpv - vecEmc).Mag() ;
188 }
189 if(distance < fRcpv + 2*delta )
92862013 190 toofar = kFALSE ;
d15a28e7 191 }
d15a28e7 192
aa0b9641 193 return distance ;
194}
195
196//____________________________________________________________________________
197TVector3 AliPHOSTrackSegmentMakerv1::PropagateToCPV(Double_t *x, Double_t *p,
198 Int_t moduleNumber) const
199{
200 // Propagate a straight-line track from the origin point x
201 // along the direction p to the CPV module moduleNumber
202 // Returns a local position of such a propagation
203
204 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
205 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
206 TVector3 moduleCenter = geom->GetCpvModuleCenter(moduleNumber);
207 TVector3 vertex(x);
208 TVector3 direction(p);
209
210 Double_t time = (moduleCenter.Mag2() - vertex.Dot(moduleCenter)) /
211 (direction.Dot(moduleCenter));
212 TVector3 globalIntersection = vertex + direction*time;
213 return geom->Global2LocalCpv(globalIntersection,moduleNumber);
d15a28e7 214}
215
7b7c1533 216//____________________________________________________________________________
217void AliPHOSTrackSegmentMakerv1::Init()
218{
219 // Make all memory allocations that are not possible in default constructor
220
88cb7938 221 AliPHOSGetter* gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
fbf811ec 222
7b7c1533 223 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
88cb7938 224 if ( !gime->TrackSegmentMaker() ) {
225 gime->PostTrackSegmentMaker(this);
226 }
7b7c1533 227}
228
8d0f3f77 229//____________________________________________________________________________
230void AliPHOSTrackSegmentMakerv1::InitParameters()
231{
e957fea8 232 //Initializes parameters
fbf811ec 233 fRcpv = 10. ;
8d0f3f77 234 fEmcFirst = 0 ;
235 fEmcLast = 0 ;
236 fCpvFirst = 0 ;
237 fCpvLast = 0 ;
238 fLinkUpArray = 0 ;
88cb7938 239 fTrackSegmentsInRun = 0 ;
eabde521 240 SetEventRange(0,-1) ;
8d0f3f77 241}
242
243
d15a28e7 244//____________________________________________________________________________
baef0810 245void AliPHOSTrackSegmentMakerv1::MakeLinks()const
d15a28e7 246{
f035f6ce 247 // Finds distances (links) between all EMC and PPSD clusters,
fbf811ec 248 // which are not further apart from each other than fRcpv
f035f6ce 249 // and sort them in accordance with this distance
9688c1dd 250
88cb7938 251 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
fbf811ec 252 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
253 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
7b7c1533 254
2731cd1e 255 fLinkUpArray->Clear() ;
2731cd1e 256
2731cd1e 257 AliPHOSRecPoint * cpv ;
92862013 258 AliPHOSEmcRecPoint * emcclu ;
28c3a259 259
d15a28e7 260 Int_t iLinkUp = 0 ;
261
28c3a259 262 Int_t iEmcRP;
2731cd1e 263 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
29b077b5 264 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
2731cd1e 265
9688c1dd 266 Bool_t toofar ;
2731cd1e 267 Int_t iCpv = 0 ;
268 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
28c3a259 269
29b077b5 270 cpv = dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(iCpv)) ;
2731cd1e 271 Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
d15a28e7 272
92862013 273 if(toofar)
88cb7938 274 break ;
fbf811ec 275 if(r < fRcpv) {
88cb7938 276 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(r, iEmcRP, iCpv) ;
28c3a259 277 }
d15a28e7 278 }
28c3a259 279 }
d15a28e7 280
9688c1dd 281 fLinkUpArray->Sort() ; //first links with smallest distances
d15a28e7 282}
28c3a259 283
d15a28e7 284//____________________________________________________________________________
2731cd1e 285void AliPHOSTrackSegmentMakerv1::MakePairs()
6ad0bfa0 286{
f035f6ce 287 // Using the previously made list of "links", we found the smallest link - i.e.
a4e98857 288 // link with the least distance between EMC and CPV and pointing to still
f035f6ce 289 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
290 // remove them from the list of "unassigned".
88cb7938 291
292 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
293
fbf811ec 294 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
295 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
88cb7938 296 TClonesArray * trackSegments = gime->TrackSegments();
9688c1dd 297
01a599c9 298 //Make arrays to mark clusters already chosen
2731cd1e 299 Int_t * emcExist = 0;
300 if(fEmcLast > fEmcFirst)
301 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
302
303 Int_t index;
304 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
305 emcExist[index] = 1 ;
306
307 Bool_t * cpvExist = 0;
308 if(fCpvLast > fCpvFirst)
309 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
310 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
311 cpvExist[index] = kTRUE ;
312
2731cd1e 313
314 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
2731cd1e 315 TIter nextUp(fLinkUpArray) ;
d15a28e7 316
d15a28e7 317 AliPHOSLink * linkUp ;
9688c1dd 318
2731cd1e 319 AliPHOSRecPoint * nullpointer = 0 ;
9688c1dd 320
29b077b5 321 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
9688c1dd 322
2731cd1e 323 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet
d15a28e7 324
2731cd1e 325 if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
88cb7938 326
327 new ((* trackSegments)[fNTrackSegments])
328 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
329 dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(linkUp->GetPpsd()))) ;
330 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
331 fNTrackSegments++ ;
332
333 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
334 //mark CPV recpoint as already used
335 cpvExist[linkUp->GetPpsd()-fCpvFirst] = kFALSE ;
7956ec10 336 } //if ppsdUp still exist
28c3a259 337 }
88cb7938 338 }
28c3a259 339
2731cd1e 340 //look through emc recPoints left without CPV/PPSD
341 if(emcExist){ //if there is emc rec point
342 Int_t iEmcRP ;
343 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
344 if(emcExist[iEmcRP] > 0 ){
88cb7938 345 new ((*trackSegments)[fNTrackSegments])
346 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
347 nullpointer) ;
348 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
349 fNTrackSegments++;
2731cd1e 350 }
d15a28e7 351 }
d15a28e7 352 }
780a31c1 353 delete [] emcExist ;
354 delete [] cpvExist ;
d15a28e7 355}
356
357//____________________________________________________________________________
eabde521 358void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
d15a28e7 359{
eabde521 360 // Steering method to perform track segment construction for events
361 // in the range from fFirstEvent to fLastEvent.
362 // This range is optionally set by SetEventRange().
363 // if fLastEvent=-1 (by default), then process events until the end.
88cb7938 364
2731cd1e 365 if(strstr(option,"tim"))
b3f97575 366 gBenchmark->Start("PHOSTSMaker");
7b7c1533 367
368 if(strstr(option,"print")) {
88cb7938 369 Print() ;
7b7c1533 370 return ;
371 }
88cb7938 372
fb43ada4 373 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
374
55ea5766 375 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
7b7c1533 376
eabde521 377 if (fLastEvent == -1) fLastEvent = gime->MaxEvent() - 1 ;
fb43ada4 378 else fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
eabde521 379 Int_t nEvents = fLastEvent - fFirstEvent + 1;
01a599c9 380
fb43ada4 381 Int_t ievent ;
382 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
aa0b9641 383 gime->Event(ievent,"R") ;
55ea5766 384 //Make some initializations
385 fNTrackSegments = 0 ;
386 fEmcFirst = 0 ;
387 fEmcLast = 0 ;
388 fCpvFirst = 0 ;
389 fCpvLast = 0 ;
88cb7938 390
391 gime->TrackSegments()->Clear();
55ea5766 392
9688c1dd 393 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
394
88cb7938 395 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
2731cd1e 396 FillOneModule() ;
2731cd1e 397 MakeLinks() ;
2731cd1e 398 MakePairs() ;
2731cd1e 399 }
28c3a259 400
90cceaf6 401 WriteTrackSegments() ;
7b7c1533 402
2731cd1e 403 if(strstr(option,"deb"))
88cb7938 404 PrintTrackSegments(option);
94de8339 405
406 //increment the total number of track segments per run
88cb7938 407 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
7b7c1533 408
2731cd1e 409 }
88cb7938 410
2731cd1e 411 if(strstr(option,"tim")){
412 gBenchmark->Stop("PHOSTSMaker");
21cd0c07 413 Info("Exec", "took %f seconds for making TS %f seconds per event",
88cb7938 414 gBenchmark->GetCpuTime("PHOSTSMaker"),
eabde521 415 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
88cb7938 416 }
417 Unload();
418}
419
420//____________________________________________________________________________
421void AliPHOSTrackSegmentMakerv1::Unload()
422{
e957fea8 423 // Unloads the task from the folder
88cb7938 424 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
425 gime->PhosLoader()->UnloadRecPoints() ;
426 gime->PhosLoader()->UnloadTracks() ;
d15a28e7 427}
7b7c1533 428
d15a28e7 429//____________________________________________________________________________
88cb7938 430void AliPHOSTrackSegmentMakerv1::Print()const
a4e98857 431{
baef0810 432 // Print TrackSegmentMaker parameters
433
21cd0c07 434 TString message("") ;
7b7c1533 435 if( strcmp(GetName(), "") != 0 ) {
21cd0c07 436 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
437 message += "Making Track segments\n" ;
438 message += "with parameters:\n" ;
439 message += " Maximal EMC - CPV (PPSD) distance (cm) %f\n" ;
440 message += "============================================\n" ;
441 Info("Print", message.Data(),fRcpv) ;
2731cd1e 442 }
443 else
21cd0c07 444 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
d15a28e7 445}
7b7c1533 446
d15a28e7 447//____________________________________________________________________________
90cceaf6 448void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
a4e98857 449{
f035f6ce 450 // Writes found TrackSegments to TreeR. Creates branches
451 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
452 // In the former branch found TrackSegments are stored, while
453 // in the latter all parameters, with which TS were made.
454 // ROOT does not allow overwriting existing branches, therefore
a4e98857 455 // first we check, if branches with the same title already exist.
f035f6ce 456 // If yes - exits without writing.
88cb7938 457
458 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
55ea5766 459
fbf811ec 460 TClonesArray * trackSegments = gime->TrackSegments() ;
55ea5766 461 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
8d0f3f77 462
88cb7938 463 TTree * treeT = gime->TreeT();
464
465 //First TS
466 Int_t bufferSize = 32000 ;
467 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
761e34c0 468 tsBranch->Fill() ;
eec3ac52 469
88cb7938 470 gime->WriteTracks("OVERWRITE");
471 gime->WriteTrackSegmentMaker("OVERWRITE");
2731cd1e 472}
98cbd830 473
98cbd830 474
2731cd1e 475//____________________________________________________________________________
a4e98857 476void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
477{
f035f6ce 478 // option deb - prints # of found TrackSegments
479 // option deb all - prints as well indexed of found RecParticles assigned to the TS
9688c1dd 480
88cb7938 481 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
21cd0c07 482
88cb7938 483 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
484 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
485 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
486
2731cd1e 487 if(strstr(option,"all")) { // printing found TS
88cb7938 488 printf("TrackSegment # EMC RP# CPV RP#\n") ;
2731cd1e 489 Int_t index;
7b7c1533 490 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
491 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
88cb7938 492 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;
2731cd1e 493 }
d15a28e7 494 }
2731cd1e 495}