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 **************************************************************************/
17 /* History of cvs commits:
20 * Revision 1.76 2005/11/17 22:29:12 hristov
21 * Faster version, no attempt to match tracks outside the PHOS acceptance
23 * Revision 1.75 2005/11/17 12:35:27 hristov
24 * Use references instead of objects. Avoid to create objects when they are not really needed
26 * Revision 1.74 2005/07/08 14:01:36 hristov
27 * Tracking in non-uniform nmagnetic field (Yu.Belikov)
29 * Revision 1.73 2005/05/28 14:19:05 schutz
30 * Compilation warnings fixed by T.P.
34 //_________________________________________________________________________
35 // Implementation version 1 of algorithm class to construct PHOS track segments
36 // Track segment for PHOS is list of
37 // EMC RecPoint + (possibly) CPV RecPoint
38 // To find TrackSegments we do the following:
39 // for each EMC RecPoints we look at
40 // CPV RecPoints in the radious fRcpv.
41 // If there is such a CPV RecPoint,
42 // we make "Link" it is just indexes of EMC and CPV RecPoint and distance
43 // between them in the PHOS plane.
44 // Then we sort "Links" and starting from the
45 // least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
47 // If there is no CPV RecPoint we make TrackSegment
48 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
49 //// In principle this class should be called from AliPHOSReconstructor, but
50 // one can use it as well in standalone mode.
52 // root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
53 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
54 // // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
55 // // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
56 // root [1] t->ExecuteTask()
57 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
58 // root [4] t->ExecuteTask("deb all time")
60 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
63 // --- ROOT system ---
65 #include "TBenchmark.h"
67 // --- Standard library ---
68 #include "Riostream.h"
69 // --- AliRoot header files ---
70 #include "AliPHOSGeometry.h"
71 #include "AliPHOSTrackSegmentMakerv1.h"
72 #include "AliPHOSTrackSegment.h"
73 #include "AliPHOSLink.h"
74 #include "AliPHOSGetter.h"
76 #include "AliESDtrack.h"
78 ClassImp( AliPHOSTrackSegmentMakerv1)
81 //____________________________________________________________________________
82 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
84 // default ctor (to be used mainly by Streamer)
87 fDefaultInit = kTRUE ;
90 //____________________________________________________________________________
91 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString alirunFileName, const TString eventFolderName)
92 :AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName)
98 fDefaultInit = kFALSE ;
102 //____________________________________________________________________________
103 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
106 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
108 delete fLinkUpArray ;
112 //____________________________________________________________________________
113 const TString AliPHOSTrackSegmentMakerv1::BranchName() const
119 //____________________________________________________________________________
120 void AliPHOSTrackSegmentMakerv1::FillOneModule()
122 // Finds first and last indexes between which
123 // clusters from one PHOS module are
125 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
127 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
128 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
131 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
132 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
133 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
137 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
139 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
140 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
145 //____________________________________________________________________________
146 Float_t AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSCpvRecPoint * cpvClu, Int_t &trackindex) const
148 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
149 // Clusters are sorted in "rows" and "columns" of width 1 cm
151 //Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
152 // if you change this value, change it as well in xxxRecPoint::Compare()
153 Float_t distance2Cpv = fRcpv ;
154 Float_t distance2Track = fRtpc ;
156 trackindex = -1 ; // closest track within fRCpv
158 TVector3 vecEmc ; // Local position of EMC recpoint
159 TVector3 vecCpv ; // Local position of CPV recpoint propagated to EMC
160 TVector3 vecDist ; // Distance between local positions of two points
162 emcClu->GetLocalPosition(vecEmc) ;
163 cpvClu->GetLocalPosition(vecCpv) ;
166 if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){
168 // Find EMC-CPV distance
169 distance2Cpv = (vecCpv - vecEmc).Mag() ;
172 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
173 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
175 Double_t rPHOS = geom->GetIPtoCrystalSurface();
177 //PH Acceptance boundaries for each PHOS module
178 Int_t nModules = geom->GetNModules();
179 Double_t * thmin = new Double_t[nModules];// theta min
180 Double_t * thmax = new Double_t[nModules];// theta max
181 Double_t * phmin = new Double_t[nModules];// phi min
182 Double_t * phmax = new Double_t[nModules];// phi max
184 for (Int_t imod=0; imod<nModules; imod++) {
185 geom->EmcModuleCoverage(imod,
186 thmin[imod],thmax[imod],
187 phmin[imod],phmax[imod]);
190 // Extrapolate the global track direction if any to CPV and find the closest track
191 Int_t nTracks = fESD->GetNumberOfTracks();
192 Int_t iClosestTrack = -1;
193 Double_t minDistance = 1e6;
194 Double_t pxyz[3], xyz[3];
195 TVector3 inPHOS; //PH Used to calculate theta and phi
199 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
200 track = fESD->GetTrack(iTrack);
202 continue ; //Skip the PHOS tracks
203 if (!track->GetXYZAt(rPHOS, fESD->GetMagneticField(), xyz))
204 continue; //track coord on the cylinder of PHOS radius
205 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
207 //PH Here one has to cut out the tracks which are not inside the PHOS
209 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
210 Double_t inPhi = inPHOS.Phi();
211 Double_t inTheta = inPHOS.Theta();
214 for (Int_t imod=0; nModules; imod++) {
215 //PH Loop on modules to check if the track enters in the acceptance
216 if (thmin[imod] < inTheta && thmax[imod] > inTheta &&
217 phmin[imod] < inPhi && phmax[imod] > inPhi) {
222 if (skip) continue; //PH Skip, if not in the PHOS acceptance
224 if (!track->GetPxPyPzAt(rPHOS, fESD->GetMagneticField(), pxyz))
225 continue; // track momentum ibid.
226 PropagateToPlane(vecDist,xyz,pxyz,"CPV",cpvClu->GetPHOSMod());
227 // Info("GetDistanceInPHOSPlane","Track %d propagation to CPV = (%f,%f,%f)",
228 // iTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
230 distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
231 // Find the closest track to the EMC recpoint
232 if (distance2Track < minDistance) {
233 minDistance = distance2Track;
234 iClosestTrack = iTrack;
243 if (iClosestTrack != -1) {
244 track = fESD->GetTrack(iClosestTrack);
245 if (track->GetPxPyPzAt(rPHOS, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
246 TVector3 vecCpvGlobal; // Global position of the CPV recpoint
247 geom->GetGlobal((AliRecPoint*)cpvClu,vecCpvGlobal);
248 for (Int_t ixyz=0; ixyz<3; ixyz++)
249 xyz[ixyz] = vecCpvGlobal[ixyz];
250 PropagateToPlane(vecDist,xyz,pxyz,"EMC",cpvClu->GetPHOSMod());
251 // Info("GetDistanceInPHOSPlane","Track %d propagation to EMC = (%f,%f,%f)",
252 // iClosestTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
254 distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
258 // // If no ESD exists, than simply find EMC-CPV distance
259 // distance = (vecCpv - vecEmc).Mag() ;
261 //if(distance2Track < fRcpv + 2*delta )
262 if(distance2Track < fRtpc )
263 trackindex = iClosestTrack ;
266 // Info("GetDistanceInPHOSPlane","cpv-emc distance is %f cm",
270 return distance2Cpv ;
273 //____________________________________________________________________________
274 void AliPHOSTrackSegmentMakerv1::PropagateToPlane(TVector3& globalIntersection,
278 Int_t moduleNumber) const
280 // Propagate a straight-line track from the origin point x
281 // along the direction p to the CPV or EMC module moduleNumber
282 // Returns a local position of such a propagation
284 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
285 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
286 TVector3 moduleCenter;
287 geom->GetModuleCenter(moduleCenter,det,moduleNumber);
288 TVector3 vertex; vertex.SetXYZ(x[0],x[1],x[2]);
289 TVector3 direction; direction.SetXYZ(p[0],p[1],p[2]);
291 // Info("PropagateToCPV","Center of the %s module %d is (%f,%f,%f)",
292 // det,moduleNumber,moduleCenter[0],moduleCenter[1],moduleCenter[2]);
294 Double_t time = (moduleCenter.Mag2() - vertex.Dot(moduleCenter)) /
295 (direction.Dot(moduleCenter));
296 vertex += direction*time;
297 geom->Global2Local(globalIntersection,vertex,moduleNumber);
300 //____________________________________________________________________________
301 void AliPHOSTrackSegmentMakerv1::Init()
303 // Make all memory allocations that are not possible in default constructor
305 AliPHOSGetter* gime = AliPHOSGetter::Instance();
307 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
309 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
310 if ( !gime->TrackSegmentMaker() ) {
311 gime->PostTrackSegmentMaker(this);
315 //____________________________________________________________________________
316 void AliPHOSTrackSegmentMakerv1::InitParameters()
318 //Initializes parameters
327 fTrackSegmentsInRun = 0 ;
328 SetEventRange(0,-1) ;
332 //____________________________________________________________________________
333 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
335 // Finds distances (links) between all EMC and CPV clusters,
336 // which are not further apart from each other than fRcpv
337 // and sort them in accordance with this distance
339 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
340 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
341 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
343 fLinkUpArray->Clear() ;
345 AliPHOSCpvRecPoint * cpv ;
346 AliPHOSEmcRecPoint * emcclu ;
351 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
352 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
356 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
358 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
360 Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, track) ;
364 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(r, iEmcRP, iCpv, track) ;
369 fLinkUpArray->Sort() ; //first links with smallest distances
372 //____________________________________________________________________________
373 void AliPHOSTrackSegmentMakerv1::MakePairs()
375 // Using the previously made list of "links", we found the smallest link - i.e.
376 // link with the least distance between EMC and CPV and pointing to still
377 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
378 // remove them from the list of "unassigned".
380 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
382 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
383 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
384 TClonesArray * trackSegments = gime->TrackSegments();
386 //Make arrays to mark clusters already chosen
387 Int_t * emcExist = 0;
388 if(fEmcLast > fEmcFirst)
389 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
392 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
393 emcExist[index] = 1 ;
395 Bool_t * cpvExist = 0;
396 if(fCpvLast > fCpvFirst)
397 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
398 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
399 cpvExist[index] = kTRUE ;
402 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
403 TIter nextUp(fLinkUpArray) ;
405 AliPHOSLink * linkUp ;
407 AliPHOSCpvRecPoint * nullpointer = 0 ;
409 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
411 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
413 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
414 new ((* trackSegments)[fNTrackSegments])
415 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
416 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
417 linkUp->GetTrack()) ;
419 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
421 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
422 //mark CPV recpoint as already used
423 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
424 } //if CpvUp still exist
428 //look through emc recPoints left without CPV
429 if(emcExist){ //if there is emc rec point
431 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
432 if(emcExist[iEmcRP] > 0 ){
433 new ((*trackSegments)[fNTrackSegments])
434 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
436 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
445 //____________________________________________________________________________
446 void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
448 // Steering method to perform track segment construction for events
449 // in the range from fFirstEvent to fLastEvent.
450 // This range is optionally set by SetEventRange().
451 // if fLastEvent=-1 (by default), then process events until the end.
453 if(strstr(option,"tim"))
454 gBenchmark->Start("PHOSTSMaker");
456 if(strstr(option,"print")) {
461 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
463 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
465 if (fLastEvent == -1)
466 fLastEvent = gime->MaxEvent() - 1 ;
468 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
469 Int_t nEvents = fLastEvent - fFirstEvent + 1;
472 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
473 gime->Event(ievent,"R") ;
474 //Make some initializations
475 fNTrackSegments = 0 ;
481 gime->TrackSegments()->Clear();
483 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
485 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
491 WriteTrackSegments() ;
493 if(strstr(option,"deb"))
494 PrintTrackSegments(option);
496 //increment the total number of track segments per run
497 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
500 if(strstr(option,"tim")){
501 gBenchmark->Stop("PHOSTSMaker");
502 Info("Exec", "took %f seconds for making TS %f seconds per event",
503 gBenchmark->GetCpuTime("PHOSTSMaker"),
504 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
506 if(fWrite) //do not unload in "on flight" mode
510 //____________________________________________________________________________
511 void AliPHOSTrackSegmentMakerv1::Unload()
513 // Unloads the task from the folder
514 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
515 gime->PhosLoader()->UnloadRecPoints() ;
516 gime->PhosLoader()->UnloadTracks() ;
519 //____________________________________________________________________________
520 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
522 // Print TrackSegmentMaker parameters
524 TString message("") ;
525 if( strcmp(GetName(), "") != 0 ) {
526 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
527 message += "Making Track segments\n" ;
528 message += "with parameters:\n" ;
529 message += " Maximal EMC - CPV distance (cm) %f\n" ;
530 message += "============================================\n" ;
531 Info("Print", message.Data(),fRcpv) ;
534 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
537 //____________________________________________________________________________
538 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
540 // Writes found TrackSegments to TreeR. Creates branches
541 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
542 // In the former branch found TrackSegments are stored, while
543 // in the latter all parameters, with which TS were made.
544 // ROOT does not allow overwriting existing branches, therefore
545 // first we check, if branches with the same title already exist.
546 // If yes - exits without writing.
548 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
550 TClonesArray * trackSegments = gime->TrackSegments() ;
551 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
553 if(fWrite){ //We write TreeT
554 TTree * treeT = gime->TreeT();
557 Int_t bufferSize = 32000 ;
558 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
561 gime->WriteTracks("OVERWRITE");
562 gime->WriteTrackSegmentMaker("OVERWRITE");
567 //____________________________________________________________________________
568 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
570 // option deb - prints # of found TrackSegments
571 // option deb all - prints as well indexed of found RecParticles assigned to the TS
573 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
575 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
576 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
577 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
579 if(strstr(option,"all")) { // printing found TS
580 printf("TrackSegment # EMC RP# CPV RP#\n") ;
582 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
583 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
584 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;