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.82 2007/03/06 06:54:48 kharlov
21 * DP:Calculation of cluster properties dep. on vertex added
23 * Revision 1.81 2007/02/05 10:02:40 kharlov
24 * Module numbering is corrected
26 * Revision 1.80 2006/08/28 10:01:56 kharlov
27 * Effective C++ warnings fixed (Timur Pocheptsov)
29 * Revision 1.79 2006/04/25 12:41:15 hristov
30 * Moving non-persistent data to AliESDfriend (Yu.Belikov)
32 * Revision 1.78 2005/11/18 13:04:51 hristov
35 * Revision 1.77 2005/11/17 23:34:36 hristov
38 * Revision 1.76 2005/11/17 22:29:12 hristov
39 * Faster version, no attempt to match tracks outside the PHOS acceptance
41 * Revision 1.75 2005/11/17 12:35:27 hristov
42 * Use references instead of objects. Avoid to create objects when they are not really needed
44 * Revision 1.74 2005/07/08 14:01:36 hristov
45 * Tracking in non-uniform nmagnetic field (Yu.Belikov)
47 * Revision 1.73 2005/05/28 14:19:05 schutz
48 * Compilation warnings fixed by T.P.
52 //_________________________________________________________________________
53 // Implementation version 1 of algorithm class to construct PHOS track segments
54 // Track segment for PHOS is list of
55 // EMC RecPoint + (possibly) CPV RecPoint
56 // To find TrackSegments we do the following:
57 // for each EMC RecPoints we look at
58 // CPV RecPoints in the radious fRcpv.
59 // If there is such a CPV RecPoint,
60 // we make "Link" it is just indexes of EMC and CPV RecPoint and distance
61 // between them in the PHOS plane.
62 // Then we sort "Links" and starting from the
63 // least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
65 // If there is no CPV RecPoint we make TrackSegment
66 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
67 //// In principle this class should be called from AliPHOSReconstructor, but
68 // one can use it as well in standalone mode.
70 // root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
71 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
72 // // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
73 // // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
74 // root [1] t->ExecuteTask()
75 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
76 // root [4] t->ExecuteTask("deb all time")
78 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
81 // --- ROOT system ---
83 #include "TBenchmark.h"
85 // --- Standard library ---
86 #include "Riostream.h"
87 // --- AliRoot header files ---
88 #include "AliPHOSGeometry.h"
89 #include "AliPHOSTrackSegmentMakerv1.h"
90 #include "AliPHOSTrackSegment.h"
91 #include "AliPHOSLink.h"
92 #include "AliPHOSGetter.h"
94 #include "AliESDtrack.h"
96 ClassImp( AliPHOSTrackSegmentMakerv1)
99 //____________________________________________________________________________
100 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() :
101 AliPHOSTrackSegmentMaker(),
113 fTrackSegmentsInRun(0)
116 // default ctor (to be used mainly by Streamer)
120 //____________________________________________________________________________
121 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString & alirunFileName, const TString & eventFolderName) :
122 AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
123 fDefaultInit(kFALSE),
134 fTrackSegmentsInRun(0)
143 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const AliPHOSTrackSegmentMakerv1 & tsm) :
144 AliPHOSTrackSegmentMaker(tsm),
145 fDefaultInit(kFALSE),
156 fTrackSegmentsInRun(0)
158 // cpy ctor: no implementation yet
159 // requested by the Coding Convention
160 Fatal("cpy ctor", "not implemented") ;
164 //____________________________________________________________________________
165 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
168 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
170 delete fLinkUpArray ;
174 //____________________________________________________________________________
175 const TString AliPHOSTrackSegmentMakerv1::BranchName() const
181 //____________________________________________________________________________
182 void AliPHOSTrackSegmentMakerv1::FillOneModule()
184 // Finds first and last indexes between which
185 // clusters from one PHOS module are
187 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
189 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
190 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
193 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
194 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
195 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
199 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
201 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
202 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
207 //____________________________________________________________________________
208 void AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
209 AliPHOSCpvRecPoint * cpvClu,
211 Float_t &dx, Float_t &dz) const
213 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
214 // Clusters are sorted in "rows" and "columns" of width 1 cm
216 // Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
217 // // if you change this value, change it as well in xxxRecPoint::Compare()
218 Float_t distance2Track = fRtpc ;
220 trackindex = -1 ; // closest track within fRCpv
222 TVector3 vecEmc ; // Local position of EMC recpoint
223 TVector3 vecP ; // Momentum direction at CPV plain
224 TVector3 vecPloc ; // Momentum direction at CPV plain
227 if(emcClu->GetPHOSMod() != cpvClu->GetPHOSMod()){
233 emcClu->GetLocalPosition(vecEmc) ;
235 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
236 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
238 Double_t xCPV,zCPV ; //EMC-projected coordinates of CPV cluster
239 TVector3 cpvGlobal; // Global position of the CPV recpoint
240 geom->GetGlobal((AliRecPoint*)cpvClu,cpvGlobal);
241 Double_t vtxCPV[3]={cpvGlobal.X(),cpvGlobal.Y(),cpvGlobal.Z()} ;
244 //if no track information available, assume straight line from IP to emcal
246 geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,xCPV,zCPV) ;
247 dx=xCPV - vecEmc.X() ;
248 dz=zCPV - vecEmc.Z() ;
252 //if there is ESD try to correct distance using TPC information on particle direct in CPV
255 Double_t rCPV = cpvGlobal.Pt() ;// Radius from IP to current point
257 // Extrapolate the global track direction if any to CPV and find the closest track
258 Int_t nTracks = fESD->GetNumberOfTracks();
259 Int_t iClosestTrack = -1;
260 Double_t minDistance = 1.e6;
266 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
267 track = fESD->GetTrack(iTrack);
268 if (!track->GetXYZAt(rCPV, fESD->GetMagneticField(), xyz))
269 continue; //track coord on the cylinder of PHOS radius
270 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
272 //Check if this track hits PHOS
273 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
274 distance2Track = inPHOS.Angle(cpvGlobal) ;
275 // Find the closest track to the CPV recpoint
276 if (distance2Track < minDistance) {
277 minDistance = distance2Track;
278 iClosestTrack = iTrack;
282 if (iClosestTrack != -1) {
283 track = fESD->GetTrack(iClosestTrack);
284 if (track->GetPxPyPzAt(rCPV, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
285 vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
287 geom->ImpactOnEmc(vtxCPV,vecP.Theta(),vecP.Phi(),dummyMod,xCPV,zCPV) ;
291 if(minDistance < fRtpc ){
292 trackindex = iClosestTrack ;
296 dx=xCPV - vecEmc.X() ;
297 dz=zCPV - vecEmc.Z() ;
305 //____________________________________________________________________________
306 void AliPHOSTrackSegmentMakerv1::Init()
308 // Make all memory allocations that are not possible in default constructor
310 AliPHOSGetter* gime = AliPHOSGetter::Instance();
312 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
314 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
315 if ( !gime->TrackSegmentMaker() ) {
316 gime->PostTrackSegmentMaker(this);
320 //____________________________________________________________________________
321 void AliPHOSTrackSegmentMakerv1::InitParameters()
323 //Initializes parameters
332 fTrackSegmentsInRun = 0 ;
333 SetEventRange(0,-1) ;
337 //____________________________________________________________________________
338 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
340 // Finds distances (links) between all EMC and CPV clusters,
341 // which are not further apart from each other than fRcpv
342 // and sort them in accordance with this distance
344 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
345 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
346 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
348 fLinkUpArray->Clear() ;
350 AliPHOSCpvRecPoint * cpv ;
351 AliPHOSEmcRecPoint * emcclu ;
356 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
357 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
361 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
363 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
366 GetDistanceInPHOSPlane(emcclu, cpv, track,dx,dz) ;
367 if(TMath::Sqrt(dx*dx+dz*dz) < fRcpv ){
368 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(dx, dz, iEmcRP, iCpv, track) ;
373 fLinkUpArray->Sort() ; //first links with smallest distances
376 //____________________________________________________________________________
377 void AliPHOSTrackSegmentMakerv1::MakePairs()
379 // Using the previously made list of "links", we found the smallest link - i.e.
380 // link with the least distance between EMC and CPV and pointing to still
381 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
382 // remove them from the list of "unassigned".
384 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
386 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
387 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
388 TClonesArray * trackSegments = gime->TrackSegments();
390 //Make arrays to mark clusters already chosen
391 Int_t * emcExist = 0;
392 if(fEmcLast > fEmcFirst)
393 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
396 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
397 emcExist[index] = 1 ;
399 Bool_t * cpvExist = 0;
400 if(fCpvLast > fCpvFirst)
401 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
402 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
403 cpvExist[index] = kTRUE ;
406 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
407 TIter nextUp(fLinkUpArray) ;
409 AliPHOSLink * linkUp ;
411 AliPHOSCpvRecPoint * nullpointer = 0 ;
413 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
415 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
417 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
419 linkUp->GetXZ(dx,dz) ;
420 new ((* trackSegments)[fNTrackSegments])
421 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
422 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
423 linkUp->GetTrack(),dx,dz) ;
425 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
427 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
428 //mark CPV recpoint as already used
429 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
430 } //if CpvUp still exist
434 //look through emc recPoints left without CPV
435 if(emcExist){ //if there is emc rec point
437 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
438 if(emcExist[iEmcRP] > 0 ){
439 new ((*trackSegments)[fNTrackSegments])
440 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
442 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
451 //____________________________________________________________________________
452 void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
454 // Steering method to perform track segment construction for events
455 // in the range from fFirstEvent to fLastEvent.
456 // This range is optionally set by SetEventRange().
457 // if fLastEvent=-1 (by default), then process events until the end.
459 if(strstr(option,"tim"))
460 gBenchmark->Start("PHOSTSMaker");
462 if(strstr(option,"print")) {
467 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
469 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
471 if (fLastEvent == -1)
472 fLastEvent = gime->MaxEvent() - 1 ;
474 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
475 Int_t nEvents = fLastEvent - fFirstEvent + 1;
478 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
479 gime->Event(ievent,"DR") ;
480 //Make some initializations
481 fNTrackSegments = 0 ;
487 gime->TrackSegments()->Clear();
492 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
494 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
500 WriteTrackSegments() ;
502 if(strstr(option,"deb"))
503 PrintTrackSegments(option);
505 //increment the total number of track segments per run
506 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
509 if(strstr(option,"tim")){
510 gBenchmark->Stop("PHOSTSMaker");
511 Info("Exec", "took %f seconds for making TS %f seconds per event",
512 gBenchmark->GetCpuTime("PHOSTSMaker"),
513 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
515 if(fWrite) //do not unload in "on flight" mode
518 //____________________________________________________________________________
519 void AliPHOSTrackSegmentMakerv1::GetVertex(void)
520 { //extract vertex either using ESD or generator
522 //Try to extract vertex from data
524 const AliESDVertex *esdVtx = fESD->GetVertex() ;
526 fVtx.SetXYZ(esdVtx->GetXv(),esdVtx->GetYv(),esdVtx->GetZv()) ;
531 AliWarning("Can not read vertex from data, use fixed \n") ;
532 fVtx.SetXYZ(0.,0.,0.) ;
535 //____________________________________________________________________________
536 void AliPHOSTrackSegmentMakerv1::EvalRecPoints(void)
537 { //calculate parameters of RecPoints using vertex and writing them
539 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
540 TClonesArray * digits = gime->Digits() ;
541 AliPHOSClusterizer * cl = gime->Clusterizer() ;
542 Double_t w0=cl->GetEmcLogWeight() ;
543 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
544 for(Int_t i=0; i<emcRecPoints->GetEntriesFast() ; i++){
545 static_cast<AliPHOSEmcRecPoint*>(emcRecPoints->At(i))->EvalAll(w0,fVtx,digits) ;
547 emcRecPoints->Sort() ;
549 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
550 Double_t w0CPV=cl->GetCpvLogWeight() ;
551 for(Int_t i=0; i<emcRecPoints->GetEntriesFast() ; i++){
552 static_cast<AliPHOSCpvRecPoint*>(cpvRecPoints->At(i))->EvalAll(w0CPV,fVtx,digits) ;
554 cpvRecPoints->Sort() ;
556 //write recaculated RecPoints
557 gime->WriteRecPoints("OVERWRITE");
558 gime->WriteClusterizer("OVERWRITE");
562 //____________________________________________________________________________
563 void AliPHOSTrackSegmentMakerv1::Unload()
565 // Unloads the task from the folder
566 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
567 gime->PhosLoader()->UnloadRecPoints() ;
568 gime->PhosLoader()->UnloadTracks() ;
571 //____________________________________________________________________________
572 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
574 // Print TrackSegmentMaker parameters
576 TString message("") ;
577 if( strcmp(GetName(), "") != 0 ) {
578 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
579 message += "Making Track segments\n" ;
580 message += "with parameters:\n" ;
581 message += " Maximal EMC - CPV distance (cm) %f\n" ;
582 message += "============================================\n" ;
583 Info("Print", message.Data(),fRcpv) ;
586 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
589 //____________________________________________________________________________
590 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
592 // Writes found TrackSegments to TreeR. Creates branches
593 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
594 // In the former branch found TrackSegments are stored, while
595 // in the latter all parameters, with which TS were made.
596 // ROOT does not allow overwriting existing branches, therefore
597 // first we check, if branches with the same title already exist.
598 // If yes - exits without writing.
600 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
602 TClonesArray * trackSegments = gime->TrackSegments() ;
603 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
605 if(fWrite){ //We write TreeT
606 TTree * treeT = gime->TreeT();
609 Int_t bufferSize = 32000 ;
610 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
613 gime->WriteTracks("OVERWRITE");
614 gime->WriteTrackSegmentMaker("OVERWRITE");
619 //____________________________________________________________________________
620 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
622 // option deb - prints # of found TrackSegments
623 // option deb all - prints as well indexed of found RecParticles assigned to the TS
625 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
627 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
628 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
629 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
631 if(strstr(option,"all")) { // printing found TS
632 printf("TrackSegment # EMC RP# CPV RP#\n") ;
634 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
635 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
636 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;