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.87 2007/05/04 14:49:29 policheh
21 * AliPHOSRecPoint inheritance from AliCluster
23 * Revision 1.86 2007/04/02 15:00:16 cvetan
24 * No more calls to gAlice in the reconstruction
26 * Revision 1.85 2007/03/28 19:18:15 kharlov
27 * RecPoints recalculation in TSM removed
29 * Revision 1.84 2007/03/07 07:01:21 hristov
30 * Fixing copy/paste erro. Additional protections
32 * Revision 1.83 2007/03/06 21:07:37 kharlov
33 * DP: xz CPV-EMC distance filled to TS
35 * Revision 1.82 2007/03/06 06:54:48 kharlov
36 * DP:Calculation of cluster properties dep. on vertex added
38 * Revision 1.81 2007/02/05 10:02:40 kharlov
39 * Module numbering is corrected
41 * Revision 1.80 2006/08/28 10:01:56 kharlov
42 * Effective C++ warnings fixed (Timur Pocheptsov)
44 * Revision 1.79 2006/04/25 12:41:15 hristov
45 * Moving non-persistent data to AliESDfriend (Yu.Belikov)
47 * Revision 1.78 2005/11/18 13:04:51 hristov
50 * Revision 1.77 2005/11/17 23:34:36 hristov
53 * Revision 1.76 2005/11/17 22:29:12 hristov
54 * Faster version, no attempt to match tracks outside the PHOS acceptance
56 * Revision 1.75 2005/11/17 12:35:27 hristov
57 * Use references instead of objects. Avoid to create objects when they are not really needed
59 * Revision 1.74 2005/07/08 14:01:36 hristov
60 * Tracking in non-uniform nmagnetic field (Yu.Belikov)
62 * Revision 1.73 2005/05/28 14:19:05 schutz
63 * Compilation warnings fixed by T.P.
67 //_________________________________________________________________________
68 // Implementation version 1 of algorithm class to construct PHOS track segments
69 // Track segment for PHOS is list of
70 // EMC RecPoint + (possibly) CPV RecPoint
71 // To find TrackSegments we do the following:
72 // for each EMC RecPoints we look at
73 // CPV RecPoints in the radious fRcpv.
74 // If there is such a CPV RecPoint,
75 // we make "Link" it is just indexes of EMC and CPV RecPoint and distance
76 // between them in the PHOS plane.
77 // Then we sort "Links" and starting from the
78 // least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
80 // If there is no CPV RecPoint we make TrackSegment
81 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
82 //// In principle this class should be called from AliPHOSReconstructor, but
83 // one can use it as well in standalone mode.
85 // root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
86 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
87 // // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
88 // // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
89 // root [1] t->ExecuteTask()
90 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
91 // root [4] t->ExecuteTask("deb all time")
93 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
96 // --- ROOT system ---
98 #include "TBenchmark.h"
100 // --- Standard library ---
101 #include "Riostream.h"
102 // --- AliRoot header files ---
103 #include "AliPHOSGeometry.h"
104 #include "AliPHOSTrackSegmentMakerv1.h"
105 #include "AliPHOSTrackSegment.h"
106 #include "AliPHOSLink.h"
107 #include "AliPHOSGetter.h"
109 #include "AliESDtrack.h"
111 ClassImp( AliPHOSTrackSegmentMakerv1)
114 //____________________________________________________________________________
115 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() :
116 AliPHOSTrackSegmentMaker(),
128 fTrackSegmentsInRun(0)
131 // default ctor (to be used mainly by Streamer)
135 //____________________________________________________________________________
136 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString & alirunFileName, const TString & eventFolderName) :
137 AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
138 fDefaultInit(kFALSE),
149 fTrackSegmentsInRun(0)
158 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const AliPHOSTrackSegmentMakerv1 & tsm) :
159 AliPHOSTrackSegmentMaker(tsm),
160 fDefaultInit(kFALSE),
171 fTrackSegmentsInRun(0)
173 // cpy ctor: no implementation yet
174 // requested by the Coding Convention
175 Fatal("cpy ctor", "not implemented") ;
179 //____________________________________________________________________________
180 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
183 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
185 delete fLinkUpArray ;
189 //____________________________________________________________________________
190 const TString AliPHOSTrackSegmentMakerv1::BranchName() const
196 //____________________________________________________________________________
197 void AliPHOSTrackSegmentMakerv1::FillOneModule()
199 // Finds first and last indexes between which
200 // clusters from one PHOS module are
202 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
204 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
205 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
208 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
209 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
210 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
214 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
216 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
217 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
222 //____________________________________________________________________________
223 void AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
224 AliPHOSCpvRecPoint * cpvClu,
226 Float_t &dx, Float_t &dz) const
228 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
229 // Clusters are sorted in "rows" and "columns" of width 1 cm
231 // Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
232 // // if you change this value, change it as well in xxxRecPoint::Compare()
233 Float_t distance2Track = fRtpc ;
235 trackindex = -1 ; // closest track within fRCpv
237 TVector3 vecEmc ; // Local position of EMC recpoint
238 TVector3 vecP ; // Momentum direction at CPV plain
239 TVector3 vecPloc ; // Momentum direction at CPV plain
242 if(emcClu->GetPHOSMod() != cpvClu->GetPHOSMod()){
248 emcClu->GetLocalPosition(vecEmc) ;
250 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
251 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
253 Double_t xCPV,zCPV ; //EMC-projected coordinates of CPV cluster
254 TVector3 cpvGlobal; // Global position of the CPV recpoint
255 geom->GetGlobalPHOS((AliPHOSRecPoint*)cpvClu,cpvGlobal);
256 Double_t vtxCPV[3]={cpvGlobal.X(),cpvGlobal.Y(),cpvGlobal.Z()} ;
260 //if no track information available, assume straight line from IP to emcal
261 geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,xCPV,zCPV) ;
262 dx=xCPV - vecEmc.X() ;
263 dz=zCPV - vecEmc.Z() ;
267 //if there is ESD try to correct distance using TPC information on particle direct in CPV
270 Double_t rCPV = cpvGlobal.Pt() ;// Radius from IP to current point
272 // Extrapolate the global track direction if any to CPV and find the closest track
273 Int_t nTracks = fESD->GetNumberOfTracks();
274 Int_t iClosestTrack = -1;
275 Double_t minDistance = 1.e6;
281 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
282 track = fESD->GetTrack(iTrack);
283 if (!track->GetXYZAt(rCPV, fESD->GetMagneticField(), xyz))
284 continue; //track coord on the cylinder of PHOS radius
285 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
287 //Check if this track hits PHOS
288 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
289 distance2Track = inPHOS.Angle(cpvGlobal) ;
290 // Find the closest track to the CPV recpoint
291 if (distance2Track < minDistance) {
292 minDistance = distance2Track;
293 iClosestTrack = iTrack;
297 if (iClosestTrack != -1) {
298 track = fESD->GetTrack(iClosestTrack);
299 if (track->GetPxPyPzAt(rCPV, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
300 vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
302 geom->ImpactOnEmc(vtxCPV,vecP.Theta(),vecP.Phi(),dummyMod,xCPV,zCPV) ;
306 if(minDistance < fRtpc ){
307 trackindex = iClosestTrack ;
311 // If the closest global track is found, calculate EMC-CPV distance from it
312 dx=xCPV - vecEmc.X() ;
313 dz=zCPV - vecEmc.Z() ;
316 // If no global track was found, just take the nearest CPV point
317 geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,xCPV,zCPV) ;
318 dx=xCPV - vecEmc.X() ;
319 dz=zCPV - vecEmc.Z() ;
323 //____________________________________________________________________________
324 void AliPHOSTrackSegmentMakerv1::Init()
326 // Make all memory allocations that are not possible in default constructor
328 AliPHOSGetter* gime = AliPHOSGetter::Instance();
330 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
332 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
333 if ( !gime->TrackSegmentMaker() ) {
334 gime->PostTrackSegmentMaker(this);
338 //____________________________________________________________________________
339 void AliPHOSTrackSegmentMakerv1::InitParameters()
341 //Initializes parameters
350 fTrackSegmentsInRun = 0 ;
351 SetEventRange(0,-1) ;
355 //____________________________________________________________________________
356 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
358 // Finds distances (links) between all EMC and CPV clusters,
359 // which are not further apart from each other than fRcpv
360 // and sort them in accordance with this distance
362 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
363 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
364 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
366 fLinkUpArray->Clear() ;
368 AliPHOSCpvRecPoint * cpv ;
369 AliPHOSEmcRecPoint * emcclu ;
374 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
375 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
379 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
381 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
384 GetDistanceInPHOSPlane(emcclu, cpv, track,dx,dz) ;
385 if(TMath::Sqrt(dx*dx+dz*dz) < fRcpv ){
386 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(dx, dz, iEmcRP, iCpv, track) ;
391 fLinkUpArray->Sort() ; //first links with smallest distances
394 //____________________________________________________________________________
395 void AliPHOSTrackSegmentMakerv1::MakePairs()
397 // Using the previously made list of "links", we found the smallest link - i.e.
398 // link with the least distance between EMC and CPV and pointing to still
399 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
400 // remove them from the list of "unassigned".
402 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
404 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
405 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
406 TClonesArray * trackSegments = gime->TrackSegments();
408 //Make arrays to mark clusters already chosen
409 Int_t * emcExist = 0;
410 if(fEmcLast > fEmcFirst)
411 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
414 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
415 emcExist[index] = 1 ;
417 Bool_t * cpvExist = 0;
418 if(fCpvLast > fCpvFirst)
419 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
420 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
421 cpvExist[index] = kTRUE ;
424 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
425 TIter nextUp(fLinkUpArray) ;
427 AliPHOSLink * linkUp ;
429 AliPHOSCpvRecPoint * nullpointer = 0 ;
431 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
433 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
435 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
437 linkUp->GetXZ(dx,dz) ;
438 new ((* trackSegments)[fNTrackSegments])
439 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
440 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
441 linkUp->GetTrack(),dx,dz) ;
443 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
445 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
446 //mark CPV recpoint as already used
447 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
448 } //if CpvUp still exist
452 //look through emc recPoints left without CPV
453 if(emcExist){ //if there is emc rec point
455 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
456 if(emcExist[iEmcRP] > 0 ){
457 new ((*trackSegments)[fNTrackSegments])
458 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
460 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
469 //____________________________________________________________________________
470 void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
472 // Steering method to perform track segment construction for events
473 // in the range from fFirstEvent to fLastEvent.
474 // This range is optionally set by SetEventRange().
475 // if fLastEvent=-1 (by default), then process events until the end.
477 if(strstr(option,"tim"))
478 gBenchmark->Start("PHOSTSMaker");
480 if(strstr(option,"print")) {
485 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
487 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
489 if (fLastEvent == -1)
490 fLastEvent = gime->MaxEvent() - 1 ;
492 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
493 Int_t nEvents = fLastEvent - fFirstEvent + 1;
496 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
497 gime->Event(ievent,"DR") ;
498 //Make some initializations
499 fNTrackSegments = 0 ;
505 gime->TrackSegments()->Clear();
507 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
509 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
515 WriteTrackSegments() ;
517 if(strstr(option,"deb"))
518 PrintTrackSegments(option);
520 //increment the total number of track segments per run
521 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
524 if(strstr(option,"tim")){
525 gBenchmark->Stop("PHOSTSMaker");
526 Info("Exec", "took %f seconds for making TS %f seconds per event",
527 gBenchmark->GetCpuTime("PHOSTSMaker"),
528 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
530 if(fWrite) //do not unload in "on flight" mode
533 //____________________________________________________________________________
534 void AliPHOSTrackSegmentMakerv1::Unload()
536 // Unloads the task from the folder
537 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
538 gime->PhosLoader()->UnloadRecPoints() ;
539 gime->PhosLoader()->UnloadTracks() ;
542 //____________________________________________________________________________
543 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
545 // Print TrackSegmentMaker parameters
547 TString message("") ;
548 if( strcmp(GetName(), "") != 0 ) {
549 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
550 message += "Making Track segments\n" ;
551 message += "with parameters:\n" ;
552 message += " Maximal EMC - CPV distance (cm) %f\n" ;
553 message += "============================================\n" ;
554 Info("Print", message.Data(),fRcpv) ;
557 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
560 //____________________________________________________________________________
561 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
563 // Writes found TrackSegments to TreeR. Creates branches
564 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
565 // In the former branch found TrackSegments are stored, while
566 // in the latter all parameters, with which TS were made.
567 // ROOT does not allow overwriting existing branches, therefore
568 // first we check, if branches with the same title already exist.
569 // If yes - exits without writing.
571 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
573 TClonesArray * trackSegments = gime->TrackSegments() ;
574 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
576 if(fWrite){ //We write TreeT
577 TTree * treeT = gime->TreeT();
580 Int_t bufferSize = 32000 ;
581 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
584 gime->WriteTracks("OVERWRITE");
585 gime->WriteTrackSegmentMaker("OVERWRITE");
590 //____________________________________________________________________________
591 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
593 // option deb - prints # of found TrackSegments
594 // option deb all - prints as well indexed of found RecParticles assigned to the TS
596 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
598 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
599 printf("nevent: %d\n", AliPHOSGetter::Instance()->EventNumber()) ;
600 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
602 if(strstr(option,"all")) { // printing found TS
603 printf("TrackSegment # EMC RP# CPV RP#\n") ;
605 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
606 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
607 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;