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