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