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.90 2007/07/11 13:43:30 hristov
21 * New class AliESDEvent, backward compatibility with the old AliESD (Christian)
23 * Revision 1.89 2007/07/03 08:13:04 kharlov
24 * Bug fix in CPV local coordinates
26 * Revision 1.88 2007/06/27 09:11:07 kharlov
27 * Bug fix for CPV-EMC distance
29 * Revision 1.87 2007/05/04 14:49:29 policheh
30 * AliPHOSRecPoint inheritance from AliCluster
32 * Revision 1.86 2007/04/02 15:00:16 cvetan
33 * No more calls to gAlice in the reconstruction
35 * Revision 1.85 2007/03/28 19:18:15 kharlov
36 * RecPoints recalculation in TSM removed
38 * Revision 1.84 2007/03/07 07:01:21 hristov
39 * Fixing copy/paste erro. Additional protections
41 * Revision 1.83 2007/03/06 21:07:37 kharlov
42 * DP: xz CPV-EMC distance filled to TS
44 * Revision 1.82 2007/03/06 06:54:48 kharlov
45 * DP:Calculation of cluster properties dep. on vertex added
47 * Revision 1.81 2007/02/05 10:02:40 kharlov
48 * Module numbering is corrected
50 * Revision 1.80 2006/08/28 10:01:56 kharlov
51 * Effective C++ warnings fixed (Timur Pocheptsov)
53 * Revision 1.79 2006/04/25 12:41:15 hristov
54 * Moving non-persistent data to AliESDfriend (Yu.Belikov)
56 * Revision 1.78 2005/11/18 13:04:51 hristov
59 * Revision 1.77 2005/11/17 23:34:36 hristov
62 * Revision 1.76 2005/11/17 22:29:12 hristov
63 * Faster version, no attempt to match tracks outside the PHOS acceptance
65 * Revision 1.75 2005/11/17 12:35:27 hristov
66 * Use references instead of objects. Avoid to create objects when they are not really needed
68 * Revision 1.74 2005/07/08 14:01:36 hristov
69 * Tracking in non-uniform nmagnetic field (Yu.Belikov)
71 * Revision 1.73 2005/05/28 14:19:05 schutz
72 * Compilation warnings fixed by T.P.
76 //_________________________________________________________________________
77 // Implementation version 1 of algorithm class to construct PHOS track segments
78 // Track segment for PHOS is list of
79 // EMC RecPoint + (possibly) CPV RecPoint
80 // To find TrackSegments we do the following:
81 // for each EMC RecPoints we look at
82 // CPV RecPoints in the radious fRcpv.
83 // If there is such a CPV RecPoint,
84 // we make "Link" it is just indexes of EMC and CPV RecPoint and distance
85 // between them in the PHOS plane.
86 // Then we sort "Links" and starting from the
87 // least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
89 // If there is no CPV RecPoint we make TrackSegment
90 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
91 //// In principle this class should be called from AliPHOSReconstructor, but
92 // one can use it as well in standalone mode.
94 // root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
95 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
96 // // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
97 // // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
98 // root [1] t->ExecuteTask()
99 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
100 // root [4] t->ExecuteTask("deb all time")
102 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
105 // --- ROOT system ---
107 #include "TBenchmark.h"
109 // --- Standard library ---
110 #include "Riostream.h"
111 // --- AliRoot header files ---
112 #include "AliPHOSGeometry.h"
113 #include "AliPHOSTrackSegmentMakerv1.h"
114 #include "AliPHOSTrackSegment.h"
115 #include "AliPHOSLink.h"
116 #include "AliPHOSGetter.h"
117 #include "AliESDEvent.h"
118 #include "AliESDtrack.h"
119 #include "AliPHOSQualAssDataMaker.h"
121 ClassImp( AliPHOSTrackSegmentMakerv1)
124 //____________________________________________________________________________
125 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() :
126 AliPHOSTrackSegmentMaker(),
139 fTrackSegmentsInRun(0)
142 // default ctor (to be used mainly by Streamer)
146 //____________________________________________________________________________
147 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString & alirunFileName, const TString & eventFolderName) :
148 AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
149 fDefaultInit(kFALSE),
161 fTrackSegmentsInRun(0)
170 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const AliPHOSTrackSegmentMakerv1 & tsm) :
171 AliPHOSTrackSegmentMaker(tsm),
172 fDefaultInit(kFALSE),
184 fTrackSegmentsInRun(0)
186 // cpy ctor: no implementation yet
187 // requested by the Coding Convention
188 Fatal("cpy ctor", "not implemented") ;
192 //____________________________________________________________________________
193 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
196 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
198 delete fLinkUpArray ;
202 //____________________________________________________________________________
203 const TString AliPHOSTrackSegmentMakerv1::BranchName() const
209 //____________________________________________________________________________
210 void AliPHOSTrackSegmentMakerv1::FillOneModule()
212 // Finds first and last indexes between which
213 // clusters from one PHOS module are
215 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
216 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
217 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
220 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
221 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
222 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
226 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
228 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
229 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
234 //____________________________________________________________________________
235 void AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
236 AliPHOSCpvRecPoint * cpvClu,
238 Float_t &dx, Float_t &dz) const
240 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
241 // Clusters are sorted in "rows" and "columns" of width 1 cm
243 // Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
244 // // if you change this value, change it as well in xxxRecPoint::Compare()
245 Float_t distance2Track = fRtpc ;
247 trackindex = -1 ; // closest track within fRCpv
249 TVector3 vecEmc ; // Local position of EMC recpoint
250 TVector3 vecP ; // Momentum direction at CPV plain
251 TVector3 vecPloc ; // Momentum direction at CPV plain
254 if(emcClu->GetPHOSMod() != cpvClu->GetPHOSMod()){
260 emcClu->GetLocalPosition(vecEmc) ;
262 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
263 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
265 Double_t xCPV,zCPV ; //EMC-projected coordinates of CPV cluster
266 TVector3 cpvGlobal; // Global position of the CPV recpoint
267 geom->GetGlobalPHOS((AliPHOSRecPoint*)cpvClu,cpvGlobal);
268 Double_t vtxCPV[3]={cpvGlobal.X(),cpvGlobal.Y(),cpvGlobal.Z()} ;
272 //if no track information available, assume straight line from IP to emcal
273 geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,zCPV,xCPV) ;
274 dx=xCPV - vecEmc.X() ;
275 dz=zCPV - vecEmc.Z() ;
279 //if there is ESD try to correct distance using TPC information on particle direct in CPV
282 Double_t rCPV = cpvGlobal.Pt() ;// Radius from IP to current point
284 // Extrapolate the global track direction if any to CPV and find the closest track
285 Int_t nTracks = fESD->GetNumberOfTracks();
286 Int_t iClosestTrack = -1;
287 Double_t minDistance = 1.e6;
293 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
294 track = fESD->GetTrack(iTrack);
295 if (!track->GetXYZAt(rCPV, fESD->GetMagneticField(), xyz))
296 continue; //track coord on the cylinder of PHOS radius
297 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
299 //Check if this track hits PHOS
300 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
301 distance2Track = inPHOS.Angle(cpvGlobal) ;
302 // Find the closest track to the CPV recpoint
303 if (distance2Track < minDistance) {
304 minDistance = distance2Track;
305 iClosestTrack = iTrack;
309 if (iClosestTrack != -1) {
310 track = fESD->GetTrack(iClosestTrack);
311 if (track->GetPxPyPzAt(rCPV, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
312 vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
314 geom->ImpactOnEmc(vtxCPV,vecP.Theta(),vecP.Phi(),dummyMod,zCPV,xCPV) ;
318 if(minDistance < fRtpc ){
319 trackindex = iClosestTrack ;
323 // If the closest global track is found, calculate EMC-CPV distance from it
324 dx=xCPV - vecEmc.X() ;
325 dz=zCPV - vecEmc.Z() ;
328 // If no global track was found, just take the nearest CPV point
329 geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,zCPV,xCPV) ;
330 dx=xCPV - vecEmc.X() ;
331 dz=zCPV - vecEmc.Z() ;
335 //____________________________________________________________________________
336 void AliPHOSTrackSegmentMakerv1::Init()
338 // Make all memory allocations that are not possible in default constructor
340 AliPHOSGetter* gime = AliPHOSGetter::Instance();
342 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
344 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
345 if ( !gime->TrackSegmentMaker() ) {
346 gime->PostTrackSegmentMaker(this);
350 //____________________________________________________________________________
351 void AliPHOSTrackSegmentMakerv1::InitParameters()
353 //Initializes parameters
362 fTrackSegmentsInRun = 0 ;
363 SetEventRange(0,-1) ;
367 //____________________________________________________________________________
368 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
370 // Finds distances (links) between all EMC and CPV clusters,
371 // which are not further apart from each other than fRcpv
372 // and sort them in accordance with this distance
374 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
375 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
376 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
378 fLinkUpArray->Clear() ;
380 AliPHOSCpvRecPoint * cpv ;
381 AliPHOSEmcRecPoint * emcclu ;
386 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
387 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
391 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
393 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
396 GetDistanceInPHOSPlane(emcclu, cpv, track,dx,dz) ;
397 if(TMath::Sqrt(dx*dx+dz*dz) < fRcpv ){
398 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(dx, dz, iEmcRP, iCpv, track) ;
403 fLinkUpArray->Sort() ; //first links with smallest distances
406 //____________________________________________________________________________
407 void AliPHOSTrackSegmentMakerv1::MakePairs()
409 // Using the previously made list of "links", we found the smallest link - i.e.
410 // link with the least distance between EMC and CPV and pointing to still
411 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
412 // remove them from the list of "unassigned".
414 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
416 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
417 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
418 TClonesArray * trackSegments = gime->TrackSegments();
420 //Make arrays to mark clusters already chosen
421 Int_t * emcExist = 0;
422 if(fEmcLast > fEmcFirst)
423 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
426 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
427 emcExist[index] = 1 ;
429 Bool_t * cpvExist = 0;
430 if(fCpvLast > fCpvFirst)
431 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
432 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
433 cpvExist[index] = kTRUE ;
436 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
437 TIter nextUp(fLinkUpArray) ;
439 AliPHOSLink * linkUp ;
441 AliPHOSCpvRecPoint * nullpointer = 0 ;
443 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
445 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
447 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
449 linkUp->GetXZ(dx,dz) ;
450 new ((* trackSegments)[fNTrackSegments])
451 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
452 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
453 linkUp->GetTrack(),dx,dz) ;
455 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
457 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
458 //mark CPV recpoint as already used
459 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
460 } //if CpvUp still exist
464 //look through emc recPoints left without CPV
465 if(emcExist){ //if there is emc rec point
467 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
468 if(emcExist[iEmcRP] > 0 ){
469 new ((*trackSegments)[fNTrackSegments])
470 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
472 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
481 //____________________________________________________________________________
482 void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
484 // Steering method to perform track segment construction for events
485 // in the range from fFirstEvent to fLastEvent.
486 // This range is optionally set by SetEventRange().
487 // if fLastEvent=-1 (by default), then process events until the end.
489 if(strstr(option,"tim"))
490 gBenchmark->Start("PHOSTSMaker");
492 if(strstr(option,"print")) {
497 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
499 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
501 if (fLastEvent == -1)
502 fLastEvent = gime->MaxEvent() - 1 ;
504 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
505 Int_t nEvents = fLastEvent - fFirstEvent + 1;
508 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
509 gime->Event(ievent,"DR") ;
510 //Make some initializations
511 fNTrackSegments = 0 ;
517 gime->TrackSegments()->Clear();
519 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
521 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
528 WriteTrackSegments() ;
530 if(strstr(option,"deb"))
531 PrintTrackSegments(option);
533 //increment the total number of track segments per run
534 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
537 if(strstr(option,"tim")){
538 gBenchmark->Stop("PHOSTSMaker");
539 Info("Exec", "took %f seconds for making TS %f seconds per event",
540 gBenchmark->GetCpuTime("PHOSTSMaker"),
541 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
543 if(fWrite) //do not unload in "on flight" mode
546 //____________________________________________________________________________
547 void AliPHOSTrackSegmentMakerv1::Unload()
549 // Unloads the task from the folder
550 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
551 gime->PhosLoader()->UnloadRecPoints() ;
552 gime->PhosLoader()->UnloadTracks() ;
555 //____________________________________________________________________________
556 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
558 // Print TrackSegmentMaker parameters
560 TString message("") ;
561 if( strcmp(GetName(), "") != 0 ) {
562 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
563 message += "Making Track segments\n" ;
564 message += "with parameters:\n" ;
565 message += " Maximal EMC - CPV distance (cm) %f\n" ;
566 message += "============================================\n" ;
567 Info("Print", message.Data(),fRcpv) ;
570 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
573 //____________________________________________________________________________
574 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
576 // Writes found TrackSegments to TreeR. Creates branches
577 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
578 // In the former branch found TrackSegments are stored, while
579 // in the latter all parameters, with which TS were made.
580 // ROOT does not allow overwriting existing branches, therefore
581 // first we check, if branches with the same title already exist.
582 // If yes - exits without writing.
584 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
586 TClonesArray * trackSegments = gime->TrackSegments() ;
587 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
589 if(fWrite){ //We write TreeT
590 TTree * treeT = gime->TreeT();
593 Int_t bufferSize = 32000 ;
594 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
597 gime->WriteTracks("OVERWRITE");
598 gime->WriteTrackSegmentMaker("OVERWRITE");
603 //____________________________________________________________________________
604 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
606 // option deb - prints # of found TrackSegments
607 // option deb all - prints as well indexed of found RecParticles assigned to the TS
609 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
611 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
612 printf("nevent: %d\n", AliPHOSGetter::Instance()->EventNumber()) ;
613 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
615 if(strstr(option,"all")) { // printing found TS
616 printf("TrackSegment # EMC RP# CPV RP#\n") ;
618 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
619 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
620 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;