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.2 2007/04/01 19:16:52 kharlov
21 * D.P.: Produce EMCTrackSegments using TPC/ITS tracks (no CPV)
26 //_________________________________________________________________________
27 // Implementation version 2 of algorithm class to construct PHOS track segments
28 // Track segment for PHOS is list of
29 // EMC RecPoint + (possibly) projection of TPC track
30 // To find TrackSegments we do the following:
31 // for each EMC RecPoints we look at
32 // TPC projections radius fRtpc.
33 // If there is such a track
34 // we make "Link" it is just indexes of EMC and TPC track and distance
35 // between them in the PHOS plane.
36 // Then we sort "Links" and starting from the
37 // least "Link" pointing to the unassined EMC and TPC assing them to
39 // If there is no TPC track we make TrackSegment
40 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
41 //// In principle this class should be called from AliPHOSReconstructor, but
42 // one can use it as well in standalone mode.
44 // root [0] AliPHOSTrackSegmentMakerv2 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
45 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
46 // // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
47 // // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
48 // root [1] t->ExecuteTask()
49 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
50 // root [4] t->ExecuteTask("deb all time")
52 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
55 // --- ROOT system ---
58 #include "TBenchmark.h"
60 // --- Standard library ---
61 #include "Riostream.h"
62 // --- AliRoot header files ---
63 #include "AliPHOSGeometry.h"
64 #include "AliPHOSTrackSegmentMakerv2.h"
65 #include "AliPHOSTrackSegment.h"
66 #include "AliPHOSLink.h"
67 #include "AliPHOSGetter.h"
69 #include "AliESDtrack.h"
71 ClassImp( AliPHOSTrackSegmentMakerv2)
74 //____________________________________________________________________________
75 AliPHOSTrackSegmentMakerv2::AliPHOSTrackSegmentMakerv2() :
76 AliPHOSTrackSegmentMaker(),
87 fTrackSegmentsInRun(0)
89 // default ctor (to be used mainly by Streamer)
93 //____________________________________________________________________________
94 AliPHOSTrackSegmentMakerv2::AliPHOSTrackSegmentMakerv2(const TString & alirunFileName, const TString & eventFolderName) :
95 AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
106 fTrackSegmentsInRun(0)
114 //____________________________________________________________________________
115 AliPHOSTrackSegmentMakerv2::AliPHOSTrackSegmentMakerv2(const AliPHOSTrackSegmentMakerv2 & tsm) :
116 AliPHOSTrackSegmentMaker(tsm),
117 fDefaultInit(kFALSE),
127 fTrackSegmentsInRun(0)
129 // cpy ctor: no implementation yet
130 // requested by the Coding Convention
131 Fatal("cpy ctor", "not implemented") ;
135 //____________________________________________________________________________
136 AliPHOSTrackSegmentMakerv2::~AliPHOSTrackSegmentMakerv2()
139 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
141 delete fLinkUpArray ;
144 //____________________________________________________________________________
145 const TString AliPHOSTrackSegmentMakerv2::BranchName() const
151 //____________________________________________________________________________
152 void AliPHOSTrackSegmentMakerv2::FillOneModule()
154 // Finds first and last indexes between which
155 // clusters from one PHOS module are
157 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
158 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
160 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
163 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
164 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
165 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
170 //Do it ones, only first time
172 Int_t nTracks = fESD->GetNumberOfTracks();
174 Int_t nPHOSmod = geom->GetNModules() ;
175 if(fTPCtracks[0].size()<(UInt_t)nTracks){
176 for(Int_t i=0; i<nPHOSmod; i++)
177 fTPCtracks[i].resize(nTracks) ;
179 for(Int_t i=0; i<5; i++)fNtpcTracks[i]=0 ;
182 //In this particular case we use fixed vertex position at zero
183 Double_t vtx[3]={0.,0.,0.} ;
186 Double_t rEMC = geom->GetIPtoCrystalSurface() ; //Use here ideal geometry
187 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
188 track = fESD->GetTrack(iTrack);
189 for(Int_t iTestMod=1; iTestMod<=nPHOSmod; iTestMod++){
190 Double_t modAngle=270.+geom->GetPHOSAngle(iTestMod) ;
191 modAngle*=TMath::Pi()/180. ;
192 track->Rotate(modAngle);
193 if (!track->GetXYZAt(rEMC, fESD->GetMagneticField(), xyz))
194 continue; //track coord on the cylinder of PHOS radius
195 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
197 //Check if this track hits PHOS
198 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
201 geom->ImpactOnEmc(vtx, inPHOS.Theta(), inPHOS.Phi(), modNum, z, x) ;
202 if(modNum==iTestMod){
203 //Mark this track as one belonging to module
204 TrackInPHOS_t &t = fTPCtracks[modNum-1][fNtpcTracks[modNum-1]++] ;
217 //____________________________________________________________________________
218 void AliPHOSTrackSegmentMakerv2::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
220 Float_t &dx, Float_t &dz) const
222 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
223 // Clusters are sorted in "rows" and "columns" of width 1 cm
225 const AliPHOSGeometry * geom = AliPHOSGeometry::GetInstance() ;
226 TVector3 emcGlobal; // Global position of the EMC recpoint
227 geom->GetGlobal((AliRecPoint*)emcClu,emcGlobal);
229 //Calculate actual distance to the PHOS surface
230 Double_t modAngle=270.+geom->GetPHOSAngle(emcClu->GetPHOSMod()) ;
231 modAngle*=TMath::Pi()/180. ;
232 Double_t rEMC = emcGlobal.X()*TMath::Cos(modAngle)+emcGlobal.Y()*TMath::Sin(modAngle) ;
233 track->Rotate(modAngle);
235 if(!track->GetXYZAt(rEMC, fESD->GetMagneticField(), xyz)){
238 return ; //track coord on the cylinder of PHOS radius
240 if((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0){
246 dx=(emcGlobal.X()-xyz[0])*TMath::Sin(modAngle)-(emcGlobal.Y()-xyz[1])*TMath::Cos(modAngle) ;
247 dx=TMath::Sign(dx,(Float_t)(emcGlobal.X()-xyz[0])) ; //set direction
248 dz=emcGlobal.Z()-xyz[2] ;
252 //____________________________________________________________________________
253 void AliPHOSTrackSegmentMakerv2::Init()
255 // Make all memory allocations that are not possible in default constructor
257 AliPHOSGetter* gime = AliPHOSGetter::Instance();
259 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
261 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
262 if ( !gime->TrackSegmentMaker() ) {
263 gime->PostTrackSegmentMaker(this);
267 //____________________________________________________________________________
268 void AliPHOSTrackSegmentMakerv2::InitParameters()
270 //Initializes parameters
276 fTrackSegmentsInRun = 0 ;
277 SetEventRange(0,-1) ;
281 //____________________________________________________________________________
282 void AliPHOSTrackSegmentMakerv2::MakeLinks()
284 // Finds distances (links) between all EMC and CPV clusters,
285 // which are not further apart from each other than fRcpv
286 // and sort them in accordance with this distance
288 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
289 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
291 fLinkUpArray->Clear() ;
293 AliPHOSEmcRecPoint * emcclu ;
298 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
299 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
301 emcclu->GetLocalPosition(vecEmc) ;
302 Int_t mod=emcclu->GetPHOSMod() ;
303 for(Int_t itr=0; itr<fNtpcTracks[mod-1] ; itr++){
304 TrackInPHOS_t &t = fTPCtracks[mod-1][itr] ;
305 //calculate raw distance
306 Double_t rawdx = t.x - vecEmc.X() ;
307 Double_t rawdz = t.z - vecEmc.Z() ;
308 if(TMath::Sqrt(rawdx*rawdx+rawdz*rawdz)<fRtpc){
310 //calculate presize difference accounting misalignment
311 GetDistanceInPHOSPlane(emcclu, t.track, dx,dz) ;
312 Int_t itrack = t.track->GetID() ;
313 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(dx, dz, iEmcRP, itrack, -1) ;
317 fLinkUpArray->Sort() ; //first links with smallest distances
320 //____________________________________________________________________________
321 void AliPHOSTrackSegmentMakerv2::MakePairs()
323 // Using the previously made list of "links", we found the smallest link - i.e.
324 // link with the least distance between EMC and CPV and pointing to still
325 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
326 // remove them from the list of "unassigned".
328 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
330 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
331 TClonesArray * trackSegments = gime->TrackSegments();
333 //Make arrays to mark clusters already chosen
334 Int_t * emcExist = 0;
335 if(fEmcLast > fEmcFirst)
336 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
339 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
340 emcExist[index] = 1 ;
343 Bool_t * tpcExist = 0;
344 Int_t nTracks = fESD->GetNumberOfTracks();
346 tpcExist = new Bool_t[nTracks] ;
348 for(index = 0; index <nTracks; index ++)
349 tpcExist[index] = 1 ;
352 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
353 TIter nextUp(fLinkUpArray) ;
355 AliPHOSLink * linkUp ;
357 AliPHOSCpvRecPoint * nullpointer = 0 ;
359 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
361 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
363 if(tpcExist[linkUp->GetCpv()]){ //Track still exist
365 linkUp->GetXZ(dx,dz) ;
366 new ((* trackSegments)[fNTrackSegments])
367 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
369 linkUp->GetTrack(),dx,dz) ;
370 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
372 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
373 //mark track as already used
374 tpcExist[linkUp->GetCpv()] = kFALSE ;
375 } //if CpvUp still exist
379 //look through emc recPoints left without CPV
380 if(emcExist){ //if there is emc rec point
382 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
383 if(emcExist[iEmcRP] > 0 ){
384 new ((*trackSegments)[fNTrackSegments])
385 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
387 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
396 //____________________________________________________________________________
397 void AliPHOSTrackSegmentMakerv2::Exec(Option_t *option)
399 // Steering method to perform track segment construction for events
400 // in the range from fFirstEvent to fLastEvent.
401 // This range is optionally set by SetEventRange().
402 // if fLastEvent=-1 (by default), then process events until the end.
404 if(strstr(option,"tim"))
405 gBenchmark->Start("PHOSTSMaker");
407 if(strstr(option,"print")) {
412 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
414 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
416 if (fLastEvent == -1)
417 fLastEvent = gime->MaxEvent() - 1 ;
419 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
420 Int_t nEvents = fLastEvent - fFirstEvent + 1;
423 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
424 gime->Event(ievent,"DR") ;
425 //Make some initializations
426 fNTrackSegments = 0 ;
430 gime->TrackSegments()->Clear();
432 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
434 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
440 WriteTrackSegments() ;
442 if(strstr(option,"deb"))
443 PrintTrackSegments(option);
445 //increment the total number of track segments per run
446 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
449 if(strstr(option,"tim")){
450 gBenchmark->Stop("PHOSTSMaker");
451 Info("Exec", "took %f seconds for making TS %f seconds per event",
452 gBenchmark->GetCpuTime("PHOSTSMaker"),
453 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
455 if(fWrite) //do not unload in "on flight" mode
458 //____________________________________________________________________________
459 void AliPHOSTrackSegmentMakerv2::Unload()
461 // Unloads the task from the folder
462 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
463 gime->PhosLoader()->UnloadRecPoints() ;
464 gime->PhosLoader()->UnloadTracks() ;
467 //____________________________________________________________________________
468 void AliPHOSTrackSegmentMakerv2::Print(const Option_t *)const
470 // Print TrackSegmentMaker parameters
472 TString message("") ;
473 if( strcmp(GetName(), "") != 0 ) {
474 message = "\n======== AliPHOSTrackSegmentMakerv2 ========\n" ;
475 message += "Making Track segments\n" ;
476 message += "with parameters:\n" ;
477 message += " Maximal EMC - TPC distance (cm) %f\n" ;
478 message += "============================================\n" ;
479 Info("Print", message.Data(),fRtpc) ;
482 Info("Print", "AliPHOSTrackSegmentMakerv2 not initialized ") ;
485 //____________________________________________________________________________
486 void AliPHOSTrackSegmentMakerv2::WriteTrackSegments()
488 // Writes found TrackSegments to TreeR. Creates branches
489 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
490 // In the former branch found TrackSegments are stored, while
491 // in the latter all parameters, with which TS were made.
492 // ROOT does not allow overwriting existing branches, therefore
493 // first we check, if branches with the same title already exist.
494 // If yes - exits without writing.
496 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
498 TClonesArray * trackSegments = gime->TrackSegments() ;
499 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
501 if(fWrite){ //We write TreeT
502 TTree * treeT = gime->TreeT();
505 Int_t bufferSize = 32000 ;
506 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
509 gime->WriteTracks("OVERWRITE");
510 gime->WriteTrackSegmentMaker("OVERWRITE");
515 //____________________________________________________________________________
516 void AliPHOSTrackSegmentMakerv2::PrintTrackSegments(Option_t * option)
518 // option deb - prints # of found TrackSegments
519 // option deb all - prints as well indexed of found RecParticles assigned to the TS
521 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
523 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
524 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
525 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
527 if(strstr(option,"all")) { // printing found TS
528 printf("TrackSegment # EMC RP# CPV RP#\n") ;
530 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
531 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
532 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetTrackIndex() ) ;