]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSTrackSegmentMakerv1.cxx
Rename AliPHOSReconstructioner to AliPHOSReconstructor
[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.
f444a19f 31//// In principle this class should be called from AliPHOSReconstructor, 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//____________________________________________________________________________
2bb500e5 129Float_t AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSCpvRecPoint * 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) ;
bfc17d18 143 cpvClu->GetLocalPosition(vecCpv) ;
2731cd1e 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;
aa0b9641 152 Double_t minDistance = 1e6;
153 Double_t pxyz[3], xyz[3];
154 AliESDtrack *track;
155 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
156 track = fESD->GetTrack(iTrack);
157 track->GetOuterXYZ(xyz); // track coord on the cylinder of PHOS radius
75408399 158 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
159 continue;
aa0b9641 160 track->GetOuterPxPyPz(pxyz); // track momentum ibid.
bfc17d18 161 vecDist = PropagateToPlane(xyz,pxyz,"CPV",cpvClu->GetPHOSMod());
162// Info("GetDistanceInPHOSPlane","Track %d propagation to CPV = (%f,%f,%f)",
163// iTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
aa0b9641 164 vecDist -= vecCpv;
165 distance = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
166 // Find the closest track to the EMC recpoint
167 if (distance < minDistance) {
168 minDistance = distance;
169 iClosestTrack = iTrack;
170 }
171 }
172
173 if (iClosestTrack != -1) {
174 track = fESD->GetTrack(iClosestTrack);
175 track->GetOuterPxPyPz(pxyz); // track momentum ibid.
176 TVector3 vecCpvGlobal; // Global position of the CPV recpoint
177 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
178 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
179 geom->GetGlobal((AliRecPoint*)cpvClu,vecCpvGlobal);
180 for (Int_t ixyz=0; ixyz<3; ixyz++)
181 xyz[ixyz] = vecCpvGlobal[ixyz];
bfc17d18 182 vecDist = PropagateToPlane(xyz,pxyz,"EMC",cpvClu->GetPHOSMod());
183// Info("GetDistanceInPHOSPlane","Track %d propagation to EMC = (%f,%f,%f)",
184// iClosestTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
aa0b9641 185 vecDist -= vecEmc;
186 distance = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
187 }
188 } else {
395f4eea 189 // If no ESD exists, than simply find EMC-CPV distance
aa0b9641 190 distance = (vecCpv - vecEmc).Mag() ;
191 }
bfc17d18 192// Info("GetDistanceInPHOSPlane","cpv-emc distance is %f cm",
193// distance);
aa0b9641 194 if(distance < fRcpv + 2*delta )
92862013 195 toofar = kFALSE ;
bfc17d18 196 }
d15a28e7 197
aa0b9641 198 return distance ;
199}
200
201//____________________________________________________________________________
bfc17d18 202TVector3 AliPHOSTrackSegmentMakerv1::PropagateToPlane(Double_t *x, Double_t *p,
203 char *det, Int_t moduleNumber) const
aa0b9641 204{
205 // Propagate a straight-line track from the origin point x
bfc17d18 206 // along the direction p to the CPV or EMC module moduleNumber
aa0b9641 207 // Returns a local position of such a propagation
208
209 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
210 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
bfc17d18 211 TVector3 moduleCenter = geom->GetModuleCenter(det,moduleNumber);
aa0b9641 212 TVector3 vertex(x);
213 TVector3 direction(p);
214
bfc17d18 215// Info("PropagateToCPV","Center of the %s module %d is (%f,%f,%f)",
216// det,moduleNumber,moduleCenter[0],moduleCenter[1],moduleCenter[2]);
395f4eea 217
aa0b9641 218 Double_t time = (moduleCenter.Mag2() - vertex.Dot(moduleCenter)) /
219 (direction.Dot(moduleCenter));
220 TVector3 globalIntersection = vertex + direction*time;
bfc17d18 221 return geom->Global2Local(globalIntersection,moduleNumber);
d15a28e7 222}
223
7b7c1533 224//____________________________________________________________________________
225void AliPHOSTrackSegmentMakerv1::Init()
226{
227 // Make all memory allocations that are not possible in default constructor
228
88cb7938 229 AliPHOSGetter* gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
fbf811ec 230
7b7c1533 231 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
88cb7938 232 if ( !gime->TrackSegmentMaker() ) {
233 gime->PostTrackSegmentMaker(this);
234 }
7b7c1533 235}
236
8d0f3f77 237//____________________________________________________________________________
238void AliPHOSTrackSegmentMakerv1::InitParameters()
239{
e957fea8 240 //Initializes parameters
fbf811ec 241 fRcpv = 10. ;
8d0f3f77 242 fEmcFirst = 0 ;
243 fEmcLast = 0 ;
244 fCpvFirst = 0 ;
245 fCpvLast = 0 ;
246 fLinkUpArray = 0 ;
88cb7938 247 fTrackSegmentsInRun = 0 ;
eabde521 248 SetEventRange(0,-1) ;
8d0f3f77 249}
250
251
d15a28e7 252//____________________________________________________________________________
baef0810 253void AliPHOSTrackSegmentMakerv1::MakeLinks()const
d15a28e7 254{
f035f6ce 255 // Finds distances (links) between all EMC and PPSD clusters,
fbf811ec 256 // which are not further apart from each other than fRcpv
f035f6ce 257 // and sort them in accordance with this distance
9688c1dd 258
88cb7938 259 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
fbf811ec 260 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
261 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
7b7c1533 262
2731cd1e 263 fLinkUpArray->Clear() ;
2731cd1e 264
2bb500e5 265 AliPHOSCpvRecPoint * cpv ;
92862013 266 AliPHOSEmcRecPoint * emcclu ;
28c3a259 267
d15a28e7 268 Int_t iLinkUp = 0 ;
269
28c3a259 270 Int_t iEmcRP;
2731cd1e 271 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
29b077b5 272 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
2731cd1e 273
9688c1dd 274 Bool_t toofar ;
2731cd1e 275 Int_t iCpv = 0 ;
276 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
28c3a259 277
2bb500e5 278 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
2731cd1e 279 Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
d15a28e7 280
92862013 281 if(toofar)
88cb7938 282 break ;
fbf811ec 283 if(r < fRcpv) {
88cb7938 284 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(r, iEmcRP, iCpv) ;
28c3a259 285 }
d15a28e7 286 }
28c3a259 287 }
d15a28e7 288
9688c1dd 289 fLinkUpArray->Sort() ; //first links with smallest distances
d15a28e7 290}
28c3a259 291
d15a28e7 292//____________________________________________________________________________
2731cd1e 293void AliPHOSTrackSegmentMakerv1::MakePairs()
6ad0bfa0 294{
f035f6ce 295 // Using the previously made list of "links", we found the smallest link - i.e.
a4e98857 296 // link with the least distance between EMC and CPV and pointing to still
f035f6ce 297 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
298 // remove them from the list of "unassigned".
88cb7938 299
300 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
301
fbf811ec 302 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
303 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
88cb7938 304 TClonesArray * trackSegments = gime->TrackSegments();
9688c1dd 305
01a599c9 306 //Make arrays to mark clusters already chosen
2731cd1e 307 Int_t * emcExist = 0;
308 if(fEmcLast > fEmcFirst)
309 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
310
311 Int_t index;
312 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
313 emcExist[index] = 1 ;
314
315 Bool_t * cpvExist = 0;
316 if(fCpvLast > fCpvFirst)
317 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
318 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
319 cpvExist[index] = kTRUE ;
320
2731cd1e 321
322 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
2731cd1e 323 TIter nextUp(fLinkUpArray) ;
d15a28e7 324
d15a28e7 325 AliPHOSLink * linkUp ;
9688c1dd 326
2bb500e5 327 AliPHOSCpvRecPoint * nullpointer = 0 ;
9688c1dd 328
29b077b5 329 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
9688c1dd 330
2731cd1e 331 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet
d15a28e7 332
2731cd1e 333 if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
88cb7938 334
335 new ((* trackSegments)[fNTrackSegments])
336 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
2bb500e5 337 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetPpsd()))) ;
88cb7938 338 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
339 fNTrackSegments++ ;
340
341 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
342 //mark CPV recpoint as already used
343 cpvExist[linkUp->GetPpsd()-fCpvFirst] = kFALSE ;
7956ec10 344 } //if ppsdUp still exist
28c3a259 345 }
88cb7938 346 }
28c3a259 347
2731cd1e 348 //look through emc recPoints left without CPV/PPSD
349 if(emcExist){ //if there is emc rec point
350 Int_t iEmcRP ;
351 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
352 if(emcExist[iEmcRP] > 0 ){
88cb7938 353 new ((*trackSegments)[fNTrackSegments])
354 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
355 nullpointer) ;
356 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
357 fNTrackSegments++;
2731cd1e 358 }
d15a28e7 359 }
d15a28e7 360 }
780a31c1 361 delete [] emcExist ;
362 delete [] cpvExist ;
d15a28e7 363}
364
365//____________________________________________________________________________
eabde521 366void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
d15a28e7 367{
eabde521 368 // Steering method to perform track segment construction for events
369 // in the range from fFirstEvent to fLastEvent.
370 // This range is optionally set by SetEventRange().
371 // if fLastEvent=-1 (by default), then process events until the end.
88cb7938 372
2731cd1e 373 if(strstr(option,"tim"))
b3f97575 374 gBenchmark->Start("PHOSTSMaker");
7b7c1533 375
376 if(strstr(option,"print")) {
88cb7938 377 Print() ;
7b7c1533 378 return ;
379 }
88cb7938 380
fb43ada4 381 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
382
55ea5766 383 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
7b7c1533 384
eabde521 385 if (fLastEvent == -1) fLastEvent = gime->MaxEvent() - 1 ;
fb43ada4 386 else fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
eabde521 387 Int_t nEvents = fLastEvent - fFirstEvent + 1;
01a599c9 388
fb43ada4 389 Int_t ievent ;
390 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
aa0b9641 391 gime->Event(ievent,"R") ;
55ea5766 392 //Make some initializations
393 fNTrackSegments = 0 ;
394 fEmcFirst = 0 ;
395 fEmcLast = 0 ;
396 fCpvFirst = 0 ;
397 fCpvLast = 0 ;
88cb7938 398
399 gime->TrackSegments()->Clear();
55ea5766 400
9688c1dd 401 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
402
88cb7938 403 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
2731cd1e 404 FillOneModule() ;
2731cd1e 405 MakeLinks() ;
2731cd1e 406 MakePairs() ;
2731cd1e 407 }
28c3a259 408
90cceaf6 409 WriteTrackSegments() ;
7b7c1533 410
2731cd1e 411 if(strstr(option,"deb"))
88cb7938 412 PrintTrackSegments(option);
94de8339 413
414 //increment the total number of track segments per run
88cb7938 415 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
7b7c1533 416
2731cd1e 417 }
88cb7938 418
2731cd1e 419 if(strstr(option,"tim")){
420 gBenchmark->Stop("PHOSTSMaker");
21cd0c07 421 Info("Exec", "took %f seconds for making TS %f seconds per event",
88cb7938 422 gBenchmark->GetCpuTime("PHOSTSMaker"),
eabde521 423 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
88cb7938 424 }
425 Unload();
426}
427
428//____________________________________________________________________________
429void AliPHOSTrackSegmentMakerv1::Unload()
430{
e957fea8 431 // Unloads the task from the folder
88cb7938 432 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
433 gime->PhosLoader()->UnloadRecPoints() ;
434 gime->PhosLoader()->UnloadTracks() ;
d15a28e7 435}
7b7c1533 436
d15a28e7 437//____________________________________________________________________________
88cb7938 438void AliPHOSTrackSegmentMakerv1::Print()const
a4e98857 439{
baef0810 440 // Print TrackSegmentMaker parameters
441
21cd0c07 442 TString message("") ;
7b7c1533 443 if( strcmp(GetName(), "") != 0 ) {
21cd0c07 444 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
445 message += "Making Track segments\n" ;
446 message += "with parameters:\n" ;
447 message += " Maximal EMC - CPV (PPSD) distance (cm) %f\n" ;
448 message += "============================================\n" ;
449 Info("Print", message.Data(),fRcpv) ;
2731cd1e 450 }
451 else
21cd0c07 452 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
d15a28e7 453}
7b7c1533 454
d15a28e7 455//____________________________________________________________________________
90cceaf6 456void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
a4e98857 457{
f035f6ce 458 // Writes found TrackSegments to TreeR. Creates branches
459 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
460 // In the former branch found TrackSegments are stored, while
461 // in the latter all parameters, with which TS were made.
462 // ROOT does not allow overwriting existing branches, therefore
a4e98857 463 // first we check, if branches with the same title already exist.
f035f6ce 464 // If yes - exits without writing.
88cb7938 465
466 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
55ea5766 467
fbf811ec 468 TClonesArray * trackSegments = gime->TrackSegments() ;
55ea5766 469 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
8d0f3f77 470
88cb7938 471 TTree * treeT = gime->TreeT();
472
473 //First TS
474 Int_t bufferSize = 32000 ;
475 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
761e34c0 476 tsBranch->Fill() ;
eec3ac52 477
88cb7938 478 gime->WriteTracks("OVERWRITE");
479 gime->WriteTrackSegmentMaker("OVERWRITE");
2731cd1e 480}
98cbd830 481
98cbd830 482
2731cd1e 483//____________________________________________________________________________
a4e98857 484void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
485{
f035f6ce 486 // option deb - prints # of found TrackSegments
487 // option deb all - prints as well indexed of found RecParticles assigned to the TS
9688c1dd 488
88cb7938 489 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
21cd0c07 490
88cb7938 491 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
492 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
493 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
494
2731cd1e 495 if(strstr(option,"all")) { // printing found TS
88cb7938 496 printf("TrackSegment # EMC RP# CPV RP#\n") ;
2731cd1e 497 Int_t index;
7b7c1533 498 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
499 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
88cb7938 500 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;
2731cd1e 501 }
d15a28e7 502 }
2731cd1e 503}