1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //_________________________________________________________________________
17 // Implementation version 1 of algorithm class to construct PHOS track segments
18 // Track segment for PHOS is list of
19 // EMC RecPoint + (possibly) CPV RecPoint
20 // To find TrackSegments we do the following:
21 // for each EMC RecPoints we look at
22 // CPV RecPoints in the radious fRcpv.
23 // If there is such a CPV RecPoint,
24 // we make "Link" it is just indexes of EMC and CPV 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
29 // If there is no CPV RecPoint we make TrackSegment
30 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
31 //// In principle this class should be called from AliPHOSReconstructor, but
32 // one can use it as well in standalone mode.
34 // root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
35 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
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")
38 // root [1] t->ExecuteTask()
39 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
40 // root [4] t->ExecuteTask("deb all time")
42 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
45 // --- ROOT system ---
47 #include "TBenchmark.h"
49 // --- Standard library ---
50 #include "Riostream.h"
51 // --- AliRoot header files ---
52 #include "AliPHOSGeometry.h"
53 #include "AliPHOSTrackSegmentMakerv1.h"
54 #include "AliPHOSTrackSegment.h"
55 #include "AliPHOSLink.h"
56 #include "AliPHOSGetter.h"
58 #include "AliESDtrack.h"
60 ClassImp( AliPHOSTrackSegmentMakerv1)
63 //____________________________________________________________________________
64 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
66 // default ctor (to be used mainly by Streamer)
69 fDefaultInit = kTRUE ;
72 //____________________________________________________________________________
73 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString alirunFileName, const TString eventFolderName)
74 :AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName)
80 fDefaultInit = kFALSE ;
84 //____________________________________________________________________________
85 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
88 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
94 //____________________________________________________________________________
95 const TString AliPHOSTrackSegmentMakerv1::BranchName() const
101 //____________________________________________________________________________
102 void AliPHOSTrackSegmentMakerv1::FillOneModule()
104 // Finds first and last indexes between which
105 // clusters from one PHOS module are
107 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
109 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
110 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
113 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
114 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
115 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
119 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
121 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
122 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
127 //____________________________________________________________________________
128 Float_t AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSCpvRecPoint * cpvClu, Int_t &trackindex) const
130 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
131 // Clusters are sorted in "rows" and "columns" of width 1 cm
133 //Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
134 // if you change this value, change it as well in xxxRecPoint::Compare()
135 Float_t distance2Cpv = fRcpv ;
136 Float_t distance2Track = fRtpc ;
138 trackindex = -1 ; // closest track within fRCpv
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
144 emcClu->GetLocalPosition(vecEmc) ;
145 cpvClu->GetLocalPosition(vecCpv) ;
148 if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){
150 // Find EMC-CPV distance
151 distance2Cpv = (vecCpv - vecEmc).Mag() ;
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;
157 Double_t minDistance = 1e6;
158 Double_t pxyz[3], xyz[3];
160 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
161 track = fESD->GetTrack(iTrack);
164 track->GetOuterXYZPHOS(xyz); // track coord on the cylinder of PHOS radius
165 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
167 track->GetOuterPxPyPzPHOS(pxyz); // track momentum ibid.
168 vecDist = PropagateToPlane(xyz,pxyz,"CPV",cpvClu->GetPHOSMod());
169 // Info("GetDistanceInPHOSPlane","Track %d propagation to CPV = (%f,%f,%f)",
170 // iTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
172 distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
173 // Find the closest track to the EMC recpoint
174 if (distance2Track < minDistance) {
175 minDistance = distance2Track;
176 iClosestTrack = iTrack;
180 if (iClosestTrack != -1) {
181 track = fESD->GetTrack(iClosestTrack);
182 track->GetOuterPxPyPzPHOS(pxyz); // track momentum ibid.
183 TVector3 vecCpvGlobal; // Global position of the CPV recpoint
184 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
185 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
186 geom->GetGlobal((AliRecPoint*)cpvClu,vecCpvGlobal);
187 for (Int_t ixyz=0; ixyz<3; ixyz++)
188 xyz[ixyz] = vecCpvGlobal[ixyz];
189 vecDist = PropagateToPlane(xyz,pxyz,"EMC",cpvClu->GetPHOSMod());
190 // Info("GetDistanceInPHOSPlane","Track %d propagation to EMC = (%f,%f,%f)",
191 // iClosestTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
193 distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
196 // // If no ESD exists, than simply find EMC-CPV distance
197 // distance = (vecCpv - vecEmc).Mag() ;
199 //if(distance2Track < fRcpv + 2*delta )
200 if(distance2Track < fRtpc )
201 trackindex = iClosestTrack ;
204 // Info("GetDistanceInPHOSPlane","cpv-emc distance is %f cm",
208 return distance2Cpv ;
211 //____________________________________________________________________________
212 TVector3 AliPHOSTrackSegmentMakerv1::PropagateToPlane(Double_t *x, Double_t *p,
213 char *det, Int_t moduleNumber) const
215 // Propagate a straight-line track from the origin point x
216 // along the direction p to the CPV or EMC module moduleNumber
217 // Returns a local position of such a propagation
219 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
220 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
221 TVector3 moduleCenter = geom->GetModuleCenter(det,moduleNumber);
223 TVector3 direction(p);
225 // Info("PropagateToCPV","Center of the %s module %d is (%f,%f,%f)",
226 // det,moduleNumber,moduleCenter[0],moduleCenter[1],moduleCenter[2]);
228 Double_t time = (moduleCenter.Mag2() - vertex.Dot(moduleCenter)) /
229 (direction.Dot(moduleCenter));
230 TVector3 globalIntersection = vertex + direction*time;
231 return geom->Global2Local(globalIntersection,moduleNumber);
234 //____________________________________________________________________________
235 void AliPHOSTrackSegmentMakerv1::Init()
237 // Make all memory allocations that are not possible in default constructor
239 AliPHOSGetter* gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
241 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
242 if ( !gime->TrackSegmentMaker() ) {
243 gime->PostTrackSegmentMaker(this);
247 //____________________________________________________________________________
248 void AliPHOSTrackSegmentMakerv1::InitParameters()
250 //Initializes parameters
258 fTrackSegmentsInRun = 0 ;
259 SetEventRange(0,-1) ;
263 //____________________________________________________________________________
264 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
266 // Finds distances (links) between all EMC and CPV clusters,
267 // which are not further apart from each other than fRcpv
268 // and sort them in accordance with this distance
270 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
271 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
272 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
274 fLinkUpArray->Clear() ;
276 AliPHOSCpvRecPoint * cpv ;
277 AliPHOSEmcRecPoint * emcclu ;
282 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
283 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
287 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
289 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
291 Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, track) ;
295 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(r, iEmcRP, iCpv, track) ;
300 fLinkUpArray->Sort() ; //first links with smallest distances
303 //____________________________________________________________________________
304 void AliPHOSTrackSegmentMakerv1::MakePairs()
306 // Using the previously made list of "links", we found the smallest link - i.e.
307 // link with the least distance between EMC and CPV and pointing to still
308 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
309 // remove them from the list of "unassigned".
311 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
313 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
314 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
315 TClonesArray * trackSegments = gime->TrackSegments();
317 //Make arrays to mark clusters already chosen
318 Int_t * emcExist = 0;
319 if(fEmcLast > fEmcFirst)
320 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
323 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
324 emcExist[index] = 1 ;
326 Bool_t * cpvExist = 0;
327 if(fCpvLast > fCpvFirst)
328 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
329 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
330 cpvExist[index] = kTRUE ;
333 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
334 TIter nextUp(fLinkUpArray) ;
336 AliPHOSLink * linkUp ;
338 AliPHOSCpvRecPoint * nullpointer = 0 ;
340 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
342 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
344 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
345 new ((* trackSegments)[fNTrackSegments])
346 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
347 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
348 linkUp->GetTrack()) ;
350 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
352 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
353 //mark CPV recpoint as already used
354 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
355 } //if CpvUp still exist
359 //look through emc recPoints left without CPV
360 if(emcExist){ //if there is emc rec point
362 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
363 if(emcExist[iEmcRP] > 0 ){
364 new ((*trackSegments)[fNTrackSegments])
365 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
367 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
376 //____________________________________________________________________________
377 void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
379 // Steering method to perform track segment construction for events
380 // in the range from fFirstEvent to fLastEvent.
381 // This range is optionally set by SetEventRange().
382 // if fLastEvent=-1 (by default), then process events until the end.
384 if(strstr(option,"tim"))
385 gBenchmark->Start("PHOSTSMaker");
387 if(strstr(option,"print")) {
392 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
394 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
396 if (fLastEvent == -1)
397 fLastEvent = gime->MaxEvent() - 1 ;
399 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
400 Int_t nEvents = fLastEvent - fFirstEvent + 1;
403 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
404 gime->Event(ievent,"R") ;
405 //Make some initializations
406 fNTrackSegments = 0 ;
412 gime->TrackSegments()->Clear();
414 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
416 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
422 WriteTrackSegments() ;
424 if(strstr(option,"deb"))
425 PrintTrackSegments(option);
427 //increment the total number of track segments per run
428 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
431 if(strstr(option,"tim")){
432 gBenchmark->Stop("PHOSTSMaker");
433 Info("Exec", "took %f seconds for making TS %f seconds per event",
434 gBenchmark->GetCpuTime("PHOSTSMaker"),
435 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
440 //____________________________________________________________________________
441 void AliPHOSTrackSegmentMakerv1::Unload()
443 // Unloads the task from the folder
444 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
445 gime->PhosLoader()->UnloadRecPoints() ;
446 gime->PhosLoader()->UnloadTracks() ;
449 //____________________________________________________________________________
450 void AliPHOSTrackSegmentMakerv1::Print()const
452 // Print TrackSegmentMaker parameters
454 TString message("") ;
455 if( strcmp(GetName(), "") != 0 ) {
456 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
457 message += "Making Track segments\n" ;
458 message += "with parameters:\n" ;
459 message += " Maximal EMC - CPV distance (cm) %f\n" ;
460 message += "============================================\n" ;
461 Info("Print", message.Data(),fRcpv) ;
464 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
467 //____________________________________________________________________________
468 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
470 // Writes found TrackSegments to TreeR. Creates branches
471 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
472 // In the former branch found TrackSegments are stored, while
473 // in the latter all parameters, with which TS were made.
474 // ROOT does not allow overwriting existing branches, therefore
475 // first we check, if branches with the same title already exist.
476 // If yes - exits without writing.
478 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
480 TClonesArray * trackSegments = gime->TrackSegments() ;
481 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
483 TTree * treeT = gime->TreeT();
486 Int_t bufferSize = 32000 ;
487 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
490 gime->WriteTracks("OVERWRITE");
491 gime->WriteTrackSegmentMaker("OVERWRITE");
495 //____________________________________________________________________________
496 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
498 // option deb - prints # of found TrackSegments
499 // option deb all - prints as well indexed of found RecParticles assigned to the TS
501 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
503 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
504 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
505 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
507 if(strstr(option,"all")) { // printing found TS
508 printf("TrackSegment # EMC RP# CPV RP#\n") ;
510 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
511 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
512 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;