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.86 2007/04/02 15:00:16 cvetan
21 * No more calls to gAlice in the reconstruction
23 * Revision 1.85 2007/03/28 19:18:15 kharlov
24 * RecPoints recalculation in TSM removed
26 * Revision 1.84 2007/03/07 07:01:21 hristov
27 * Fixing copy/paste erro. Additional protections
29 * Revision 1.83 2007/03/06 21:07:37 kharlov
30 * DP: xz CPV-EMC distance filled to TS
32 * Revision 1.82 2007/03/06 06:54:48 kharlov
33 * DP:Calculation of cluster properties dep. on vertex added
35 * Revision 1.81 2007/02/05 10:02:40 kharlov
36 * Module numbering is corrected
38 * Revision 1.80 2006/08/28 10:01:56 kharlov
39 * Effective C++ warnings fixed (Timur Pocheptsov)
41 * Revision 1.79 2006/04/25 12:41:15 hristov
42 * Moving non-persistent data to AliESDfriend (Yu.Belikov)
44 * Revision 1.78 2005/11/18 13:04:51 hristov
47 * Revision 1.77 2005/11/17 23:34:36 hristov
50 * Revision 1.76 2005/11/17 22:29:12 hristov
51 * Faster version, no attempt to match tracks outside the PHOS acceptance
53 * Revision 1.75 2005/11/17 12:35:27 hristov
54 * Use references instead of objects. Avoid to create objects when they are not really needed
56 * Revision 1.74 2005/07/08 14:01:36 hristov
57 * Tracking in non-uniform nmagnetic field (Yu.Belikov)
59 * Revision 1.73 2005/05/28 14:19:05 schutz
60 * Compilation warnings fixed by T.P.
64 //_________________________________________________________________________
65 // Implementation version 1 of algorithm class to construct PHOS track segments
66 // Track segment for PHOS is list of
67 // EMC RecPoint + (possibly) CPV RecPoint
68 // To find TrackSegments we do the following:
69 // for each EMC RecPoints we look at
70 // CPV RecPoints in the radious fRcpv.
71 // If there is such a CPV RecPoint,
72 // we make "Link" it is just indexes of EMC and CPV RecPoint and distance
73 // between them in the PHOS plane.
74 // Then we sort "Links" and starting from the
75 // least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
77 // If there is no CPV RecPoint we make TrackSegment
78 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
79 //// In principle this class should be called from AliPHOSReconstructor, but
80 // one can use it as well in standalone mode.
82 // root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
83 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
84 // // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
85 // // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
86 // root [1] t->ExecuteTask()
87 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
88 // root [4] t->ExecuteTask("deb all time")
90 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
93 // --- ROOT system ---
95 #include "TBenchmark.h"
97 // --- Standard library ---
98 #include "Riostream.h"
99 // --- AliRoot header files ---
100 #include "AliPHOSGeometry.h"
101 #include "AliPHOSTrackSegmentMakerv1.h"
102 #include "AliPHOSTrackSegment.h"
103 #include "AliPHOSLink.h"
104 #include "AliPHOSGetter.h"
106 #include "AliESDtrack.h"
108 ClassImp( AliPHOSTrackSegmentMakerv1)
111 //____________________________________________________________________________
112 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() :
113 AliPHOSTrackSegmentMaker(),
125 fTrackSegmentsInRun(0)
128 // default ctor (to be used mainly by Streamer)
132 //____________________________________________________________________________
133 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString & alirunFileName, const TString & eventFolderName) :
134 AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
135 fDefaultInit(kFALSE),
146 fTrackSegmentsInRun(0)
155 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const AliPHOSTrackSegmentMakerv1 & tsm) :
156 AliPHOSTrackSegmentMaker(tsm),
157 fDefaultInit(kFALSE),
168 fTrackSegmentsInRun(0)
170 // cpy ctor: no implementation yet
171 // requested by the Coding Convention
172 Fatal("cpy ctor", "not implemented") ;
176 //____________________________________________________________________________
177 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
180 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
182 delete fLinkUpArray ;
186 //____________________________________________________________________________
187 const TString AliPHOSTrackSegmentMakerv1::BranchName() const
193 //____________________________________________________________________________
194 void AliPHOSTrackSegmentMakerv1::FillOneModule()
196 // Finds first and last indexes between which
197 // clusters from one PHOS module are
199 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
201 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
202 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
205 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
206 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
207 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
211 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
213 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
214 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
219 //____________________________________________________________________________
220 void AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
221 AliPHOSCpvRecPoint * cpvClu,
223 Float_t &dx, Float_t &dz) const
225 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
226 // Clusters are sorted in "rows" and "columns" of width 1 cm
228 // Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
229 // // if you change this value, change it as well in xxxRecPoint::Compare()
230 Float_t distance2Track = fRtpc ;
232 trackindex = -1 ; // closest track within fRCpv
234 TVector3 vecEmc ; // Local position of EMC recpoint
235 TVector3 vecP ; // Momentum direction at CPV plain
236 TVector3 vecPloc ; // Momentum direction at CPV plain
239 if(emcClu->GetPHOSMod() != cpvClu->GetPHOSMod()){
245 emcClu->GetLocalPosition(vecEmc) ;
247 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
248 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
250 Double_t xCPV,zCPV ; //EMC-projected coordinates of CPV cluster
251 TVector3 cpvGlobal; // Global position of the CPV recpoint
252 geom->GetGlobalPHOS((AliPHOSRecPoint*)cpvClu,cpvGlobal);
253 Double_t vtxCPV[3]={cpvGlobal.X(),cpvGlobal.Y(),cpvGlobal.Z()} ;
256 //if no track information available, assume straight line from IP to emcal
258 geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,xCPV,zCPV) ;
259 dx=xCPV - vecEmc.X() ;
260 dz=zCPV - vecEmc.Z() ;
264 //if there is ESD try to correct distance using TPC information on particle direct in CPV
267 Double_t rCPV = cpvGlobal.Pt() ;// Radius from IP to current point
269 // Extrapolate the global track direction if any to CPV and find the closest track
270 Int_t nTracks = fESD->GetNumberOfTracks();
271 Int_t iClosestTrack = -1;
272 Double_t minDistance = 1.e6;
278 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
279 track = fESD->GetTrack(iTrack);
280 if (!track->GetXYZAt(rCPV, fESD->GetMagneticField(), xyz))
281 continue; //track coord on the cylinder of PHOS radius
282 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
284 //Check if this track hits PHOS
285 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
286 distance2Track = inPHOS.Angle(cpvGlobal) ;
287 // Find the closest track to the CPV recpoint
288 if (distance2Track < minDistance) {
289 minDistance = distance2Track;
290 iClosestTrack = iTrack;
294 if (iClosestTrack != -1) {
295 track = fESD->GetTrack(iClosestTrack);
296 if (track->GetPxPyPzAt(rCPV, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
297 vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
299 geom->ImpactOnEmc(vtxCPV,vecP.Theta(),vecP.Phi(),dummyMod,xCPV,zCPV) ;
303 if(minDistance < fRtpc ){
304 trackindex = iClosestTrack ;
308 dx=xCPV - vecEmc.X() ;
309 dz=zCPV - vecEmc.Z() ;
317 //____________________________________________________________________________
318 void AliPHOSTrackSegmentMakerv1::Init()
320 // Make all memory allocations that are not possible in default constructor
322 AliPHOSGetter* gime = AliPHOSGetter::Instance();
324 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
326 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
327 if ( !gime->TrackSegmentMaker() ) {
328 gime->PostTrackSegmentMaker(this);
332 //____________________________________________________________________________
333 void AliPHOSTrackSegmentMakerv1::InitParameters()
335 //Initializes parameters
344 fTrackSegmentsInRun = 0 ;
345 SetEventRange(0,-1) ;
349 //____________________________________________________________________________
350 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
352 // Finds distances (links) between all EMC and CPV clusters,
353 // which are not further apart from each other than fRcpv
354 // and sort them in accordance with this distance
356 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
357 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
358 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
360 fLinkUpArray->Clear() ;
362 AliPHOSCpvRecPoint * cpv ;
363 AliPHOSEmcRecPoint * emcclu ;
368 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
369 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
373 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
375 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
378 GetDistanceInPHOSPlane(emcclu, cpv, track,dx,dz) ;
379 if(TMath::Sqrt(dx*dx+dz*dz) < fRcpv ){
380 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(dx, dz, iEmcRP, iCpv, track) ;
385 fLinkUpArray->Sort() ; //first links with smallest distances
388 //____________________________________________________________________________
389 void AliPHOSTrackSegmentMakerv1::MakePairs()
391 // Using the previously made list of "links", we found the smallest link - i.e.
392 // link with the least distance between EMC and CPV and pointing to still
393 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
394 // remove them from the list of "unassigned".
396 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
398 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
399 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
400 TClonesArray * trackSegments = gime->TrackSegments();
402 //Make arrays to mark clusters already chosen
403 Int_t * emcExist = 0;
404 if(fEmcLast > fEmcFirst)
405 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
408 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
409 emcExist[index] = 1 ;
411 Bool_t * cpvExist = 0;
412 if(fCpvLast > fCpvFirst)
413 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
414 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
415 cpvExist[index] = kTRUE ;
418 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
419 TIter nextUp(fLinkUpArray) ;
421 AliPHOSLink * linkUp ;
423 AliPHOSCpvRecPoint * nullpointer = 0 ;
425 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
427 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
429 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
431 linkUp->GetXZ(dx,dz) ;
432 new ((* trackSegments)[fNTrackSegments])
433 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
434 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
435 linkUp->GetTrack(),dx,dz) ;
437 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
439 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
440 //mark CPV recpoint as already used
441 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
442 } //if CpvUp still exist
446 //look through emc recPoints left without CPV
447 if(emcExist){ //if there is emc rec point
449 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
450 if(emcExist[iEmcRP] > 0 ){
451 new ((*trackSegments)[fNTrackSegments])
452 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
454 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
463 //____________________________________________________________________________
464 void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
466 // Steering method to perform track segment construction for events
467 // in the range from fFirstEvent to fLastEvent.
468 // This range is optionally set by SetEventRange().
469 // if fLastEvent=-1 (by default), then process events until the end.
471 if(strstr(option,"tim"))
472 gBenchmark->Start("PHOSTSMaker");
474 if(strstr(option,"print")) {
479 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
481 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
483 if (fLastEvent == -1)
484 fLastEvent = gime->MaxEvent() - 1 ;
486 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
487 Int_t nEvents = fLastEvent - fFirstEvent + 1;
490 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
491 gime->Event(ievent,"DR") ;
492 //Make some initializations
493 fNTrackSegments = 0 ;
499 gime->TrackSegments()->Clear();
501 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
503 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
509 WriteTrackSegments() ;
511 if(strstr(option,"deb"))
512 PrintTrackSegments(option);
514 //increment the total number of track segments per run
515 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
518 if(strstr(option,"tim")){
519 gBenchmark->Stop("PHOSTSMaker");
520 Info("Exec", "took %f seconds for making TS %f seconds per event",
521 gBenchmark->GetCpuTime("PHOSTSMaker"),
522 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
524 if(fWrite) //do not unload in "on flight" mode
527 //____________________________________________________________________________
528 void AliPHOSTrackSegmentMakerv1::Unload()
530 // Unloads the task from the folder
531 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
532 gime->PhosLoader()->UnloadRecPoints() ;
533 gime->PhosLoader()->UnloadTracks() ;
536 //____________________________________________________________________________
537 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
539 // Print TrackSegmentMaker parameters
541 TString message("") ;
542 if( strcmp(GetName(), "") != 0 ) {
543 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
544 message += "Making Track segments\n" ;
545 message += "with parameters:\n" ;
546 message += " Maximal EMC - CPV distance (cm) %f\n" ;
547 message += "============================================\n" ;
548 Info("Print", message.Data(),fRcpv) ;
551 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
554 //____________________________________________________________________________
555 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
557 // Writes found TrackSegments to TreeR. Creates branches
558 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
559 // In the former branch found TrackSegments are stored, while
560 // in the latter all parameters, with which TS were made.
561 // ROOT does not allow overwriting existing branches, therefore
562 // first we check, if branches with the same title already exist.
563 // If yes - exits without writing.
565 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
567 TClonesArray * trackSegments = gime->TrackSegments() ;
568 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
570 if(fWrite){ //We write TreeT
571 TTree * treeT = gime->TreeT();
574 Int_t bufferSize = 32000 ;
575 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
578 gime->WriteTracks("OVERWRITE");
579 gime->WriteTrackSegmentMaker("OVERWRITE");
584 //____________________________________________________________________________
585 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
587 // option deb - prints # of found TrackSegments
588 // option deb all - prints as well indexed of found RecParticles assigned to the TS
590 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
592 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
593 printf("nevent: %d\n", AliPHOSGetter::Instance()->EventNumber()) ;
594 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
596 if(strstr(option,"all")) { // printing found TS
597 printf("TrackSegment # EMC RP# CPV RP#\n") ;
599 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
600 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
601 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;