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.83 2007/03/06 21:07:37 kharlov
21 * DP: xz CPV-EMC distance filled to TS
23 * Revision 1.82 2007/03/06 06:54:48 kharlov
24 * DP:Calculation of cluster properties dep. on vertex added
26 * Revision 1.81 2007/02/05 10:02:40 kharlov
27 * Module numbering is corrected
29 * Revision 1.80 2006/08/28 10:01:56 kharlov
30 * Effective C++ warnings fixed (Timur Pocheptsov)
32 * Revision 1.79 2006/04/25 12:41:15 hristov
33 * Moving non-persistent data to AliESDfriend (Yu.Belikov)
35 * Revision 1.78 2005/11/18 13:04:51 hristov
38 * Revision 1.77 2005/11/17 23:34:36 hristov
41 * Revision 1.76 2005/11/17 22:29:12 hristov
42 * Faster version, no attempt to match tracks outside the PHOS acceptance
44 * Revision 1.75 2005/11/17 12:35:27 hristov
45 * Use references instead of objects. Avoid to create objects when they are not really needed
47 * Revision 1.74 2005/07/08 14:01:36 hristov
48 * Tracking in non-uniform nmagnetic field (Yu.Belikov)
50 * Revision 1.73 2005/05/28 14:19:05 schutz
51 * Compilation warnings fixed by T.P.
55 //_________________________________________________________________________
56 // Implementation version 1 of algorithm class to construct PHOS track segments
57 // Track segment for PHOS is list of
58 // EMC RecPoint + (possibly) CPV RecPoint
59 // To find TrackSegments we do the following:
60 // for each EMC RecPoints we look at
61 // CPV RecPoints in the radious fRcpv.
62 // If there is such a CPV RecPoint,
63 // we make "Link" it is just indexes of EMC and CPV RecPoint and distance
64 // between them in the PHOS plane.
65 // Then we sort "Links" and starting from the
66 // least "Link" pointing to the unassined EMC and CPV RecPoints assing them to
68 // If there is no CPV RecPoint we make TrackSegment
69 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
70 //// In principle this class should be called from AliPHOSReconstructor, but
71 // one can use it as well in standalone mode.
73 // root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
74 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
75 // // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
76 // // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")
77 // root [1] t->ExecuteTask()
78 // root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
79 // root [4] t->ExecuteTask("deb all time")
81 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH)
84 // --- ROOT system ---
86 #include "TBenchmark.h"
88 // --- Standard library ---
89 #include "Riostream.h"
90 // --- AliRoot header files ---
91 #include "AliPHOSGeometry.h"
92 #include "AliPHOSTrackSegmentMakerv1.h"
93 #include "AliPHOSTrackSegment.h"
94 #include "AliPHOSLink.h"
95 #include "AliPHOSGetter.h"
97 #include "AliESDtrack.h"
99 ClassImp( AliPHOSTrackSegmentMakerv1)
102 //____________________________________________________________________________
103 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() :
104 AliPHOSTrackSegmentMaker(),
116 fTrackSegmentsInRun(0)
119 // default ctor (to be used mainly by Streamer)
123 //____________________________________________________________________________
124 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString & alirunFileName, const TString & eventFolderName) :
125 AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
126 fDefaultInit(kFALSE),
137 fTrackSegmentsInRun(0)
146 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const AliPHOSTrackSegmentMakerv1 & tsm) :
147 AliPHOSTrackSegmentMaker(tsm),
148 fDefaultInit(kFALSE),
159 fTrackSegmentsInRun(0)
161 // cpy ctor: no implementation yet
162 // requested by the Coding Convention
163 Fatal("cpy ctor", "not implemented") ;
167 //____________________________________________________________________________
168 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
171 // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
173 delete fLinkUpArray ;
177 //____________________________________________________________________________
178 const TString AliPHOSTrackSegmentMakerv1::BranchName() const
184 //____________________________________________________________________________
185 void AliPHOSTrackSegmentMakerv1::FillOneModule()
187 // Finds first and last indexes between which
188 // clusters from one PHOS module are
190 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
192 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
193 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
196 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
197 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
198 ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule );
202 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
204 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
205 ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule );
210 //____________________________________________________________________________
211 void AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
212 AliPHOSCpvRecPoint * cpvClu,
214 Float_t &dx, Float_t &dz) const
216 // Calculates the distance between the EMC RecPoint and the CPV RecPoint
217 // Clusters are sorted in "rows" and "columns" of width 1 cm
219 // Float_t delta = 1 ; // Width of the rows in sorting of RecPoints (in cm)
220 // // if you change this value, change it as well in xxxRecPoint::Compare()
221 Float_t distance2Track = fRtpc ;
223 trackindex = -1 ; // closest track within fRCpv
225 TVector3 vecEmc ; // Local position of EMC recpoint
226 TVector3 vecP ; // Momentum direction at CPV plain
227 TVector3 vecPloc ; // Momentum direction at CPV plain
230 if(emcClu->GetPHOSMod() != cpvClu->GetPHOSMod()){
236 emcClu->GetLocalPosition(vecEmc) ;
238 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
239 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
241 Double_t xCPV,zCPV ; //EMC-projected coordinates of CPV cluster
242 TVector3 cpvGlobal; // Global position of the CPV recpoint
243 geom->GetGlobal((AliRecPoint*)cpvClu,cpvGlobal);
244 Double_t vtxCPV[3]={cpvGlobal.X(),cpvGlobal.Y(),cpvGlobal.Z()} ;
247 //if no track information available, assume straight line from IP to emcal
249 geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,xCPV,zCPV) ;
250 dx=xCPV - vecEmc.X() ;
251 dz=zCPV - vecEmc.Z() ;
255 //if there is ESD try to correct distance using TPC information on particle direct in CPV
258 Double_t rCPV = cpvGlobal.Pt() ;// Radius from IP to current point
260 // Extrapolate the global track direction if any to CPV and find the closest track
261 Int_t nTracks = fESD->GetNumberOfTracks();
262 Int_t iClosestTrack = -1;
263 Double_t minDistance = 1.e6;
269 for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
270 track = fESD->GetTrack(iTrack);
271 if (!track->GetXYZAt(rCPV, fESD->GetMagneticField(), xyz))
272 continue; //track coord on the cylinder of PHOS radius
273 if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
275 //Check if this track hits PHOS
276 inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
277 distance2Track = inPHOS.Angle(cpvGlobal) ;
278 // Find the closest track to the CPV recpoint
279 if (distance2Track < minDistance) {
280 minDistance = distance2Track;
281 iClosestTrack = iTrack;
285 if (iClosestTrack != -1) {
286 track = fESD->GetTrack(iClosestTrack);
287 if (track->GetPxPyPzAt(rCPV, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
288 vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
290 geom->ImpactOnEmc(vtxCPV,vecP.Theta(),vecP.Phi(),dummyMod,xCPV,zCPV) ;
294 if(minDistance < fRtpc ){
295 trackindex = iClosestTrack ;
299 dx=xCPV - vecEmc.X() ;
300 dz=zCPV - vecEmc.Z() ;
308 //____________________________________________________________________________
309 void AliPHOSTrackSegmentMakerv1::Init()
311 // Make all memory allocations that are not possible in default constructor
313 AliPHOSGetter* gime = AliPHOSGetter::Instance();
315 gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
317 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
318 if ( !gime->TrackSegmentMaker() ) {
319 gime->PostTrackSegmentMaker(this);
323 //____________________________________________________________________________
324 void AliPHOSTrackSegmentMakerv1::InitParameters()
326 //Initializes parameters
335 fTrackSegmentsInRun = 0 ;
336 SetEventRange(0,-1) ;
340 //____________________________________________________________________________
341 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
343 // Finds distances (links) between all EMC and CPV clusters,
344 // which are not further apart from each other than fRcpv
345 // and sort them in accordance with this distance
347 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
348 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
349 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
351 fLinkUpArray->Clear() ;
353 AliPHOSCpvRecPoint * cpv ;
354 AliPHOSEmcRecPoint * emcclu ;
359 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
360 emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
364 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
366 cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
369 GetDistanceInPHOSPlane(emcclu, cpv, track,dx,dz) ;
370 if(TMath::Sqrt(dx*dx+dz*dz) < fRcpv ){
371 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(dx, dz, iEmcRP, iCpv, track) ;
376 fLinkUpArray->Sort() ; //first links with smallest distances
379 //____________________________________________________________________________
380 void AliPHOSTrackSegmentMakerv1::MakePairs()
382 // Using the previously made list of "links", we found the smallest link - i.e.
383 // link with the least distance between EMC and CPV and pointing to still
384 // unassigned RecParticles. We assign these RecPoints to TrackSegment and
385 // remove them from the list of "unassigned".
387 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
389 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
390 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
391 TClonesArray * trackSegments = gime->TrackSegments();
393 //Make arrays to mark clusters already chosen
394 Int_t * emcExist = 0;
395 if(fEmcLast > fEmcFirst)
396 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
399 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
400 emcExist[index] = 1 ;
402 Bool_t * cpvExist = 0;
403 if(fCpvLast > fCpvFirst)
404 cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
405 for(index = 0; index <fCpvLast-fCpvFirst; index ++)
406 cpvExist[index] = kTRUE ;
409 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
410 TIter nextUp(fLinkUpArray) ;
412 AliPHOSLink * linkUp ;
414 AliPHOSCpvRecPoint * nullpointer = 0 ;
416 while ( (linkUp = static_cast<AliPHOSLink *>(nextUp()) ) ){
418 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
420 if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
422 linkUp->GetXZ(dx,dz) ;
423 new ((* trackSegments)[fNTrackSegments])
424 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) ,
425 dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) ,
426 linkUp->GetTrack(),dx,dz) ;
428 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
430 emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc that Cpv was found
431 //mark CPV recpoint as already used
432 cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
433 } //if CpvUp still exist
437 //look through emc recPoints left without CPV
438 if(emcExist){ //if there is emc rec point
440 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
441 if(emcExist[iEmcRP] > 0 ){
442 new ((*trackSegments)[fNTrackSegments])
443 AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)),
445 (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
454 //____________________________________________________________________________
455 void AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
457 // Steering method to perform track segment construction for events
458 // in the range from fFirstEvent to fLastEvent.
459 // This range is optionally set by SetEventRange().
460 // if fLastEvent=-1 (by default), then process events until the end.
462 if(strstr(option,"tim"))
463 gBenchmark->Start("PHOSTSMaker");
465 if(strstr(option,"print")) {
470 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
472 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
474 if (fLastEvent == -1)
475 fLastEvent = gime->MaxEvent() - 1 ;
477 fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
478 Int_t nEvents = fLastEvent - fFirstEvent + 1;
481 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
482 gime->Event(ievent,"DR") ;
483 //Make some initializations
484 fNTrackSegments = 0 ;
490 gime->TrackSegments()->Clear();
495 // if(!ReadRecPoints(ievent)) continue; //reads RecPoints for event ievent
497 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
503 WriteTrackSegments() ;
505 if(strstr(option,"deb"))
506 PrintTrackSegments(option);
508 //increment the total number of track segments per run
509 fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ;
512 if(strstr(option,"tim")){
513 gBenchmark->Stop("PHOSTSMaker");
514 Info("Exec", "took %f seconds for making TS %f seconds per event",
515 gBenchmark->GetCpuTime("PHOSTSMaker"),
516 gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
518 if(fWrite) //do not unload in "on flight" mode
521 //____________________________________________________________________________
522 void AliPHOSTrackSegmentMakerv1::GetVertex(void)
523 { //extract vertex either using ESD or generator
525 //Try to extract vertex from data
527 const AliESDVertex *esdVtx = fESD->GetVertex() ;
529 fVtx.SetXYZ(esdVtx->GetXv(),esdVtx->GetYv(),esdVtx->GetZv()) ;
534 AliWarning("Can not read vertex from data, use fixed \n") ;
535 fVtx.SetXYZ(0.,0.,0.) ;
538 //____________________________________________________________________________
539 void AliPHOSTrackSegmentMakerv1::EvalRecPoints(void)
540 { //calculate parameters of RecPoints using vertex and writing them
542 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
543 TClonesArray * digits = gime->Digits() ;
544 AliPHOSClusterizer * cl = gime->Clusterizer() ;
545 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
547 AliError("No CPV rec. points!");
550 Double_t w0=cl->GetEmcLogWeight() ;
551 for(Int_t i=0; i<emcRecPoints->GetEntriesFast() ; i++){
552 AliPHOSEmcRecPoint * point = static_cast<AliPHOSEmcRecPoint*>(emcRecPoints->At(i));
553 if (point) point->EvalAll(w0,fVtx,digits) ;
555 AliError(Form("No AliPHOSEmcRecPoint is found at %d",i));
558 emcRecPoints->Sort() ;
561 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
563 AliError("No CPV rec. points!");
566 Double_t w0CPV=cl->GetCpvLogWeight() ;
567 for(Int_t i=0; i<cpvRecPoints->GetEntriesFast() ; i++){
568 AliPHOSCpvRecPoint * point = static_cast<AliPHOSCpvRecPoint*>(cpvRecPoints->At(i));
569 if (point) point->EvalAll(w0CPV,fVtx,digits) ;
571 AliError(Form("No AliPHOSCpvRecPoint is found at %d",i));
574 cpvRecPoints->Sort() ;
577 //write recaculated RecPoints
578 gime->WriteRecPoints("OVERWRITE");
579 gime->WriteClusterizer("OVERWRITE");
583 //____________________________________________________________________________
584 void AliPHOSTrackSegmentMakerv1::Unload()
586 // Unloads the task from the folder
587 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
588 gime->PhosLoader()->UnloadRecPoints() ;
589 gime->PhosLoader()->UnloadTracks() ;
592 //____________________________________________________________________________
593 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
595 // Print TrackSegmentMaker parameters
597 TString message("") ;
598 if( strcmp(GetName(), "") != 0 ) {
599 message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ;
600 message += "Making Track segments\n" ;
601 message += "with parameters:\n" ;
602 message += " Maximal EMC - CPV distance (cm) %f\n" ;
603 message += "============================================\n" ;
604 Info("Print", message.Data(),fRcpv) ;
607 Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
610 //____________________________________________________________________________
611 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
613 // Writes found TrackSegments to TreeR. Creates branches
614 // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
615 // In the former branch found TrackSegments are stored, while
616 // in the latter all parameters, with which TS were made.
617 // ROOT does not allow overwriting existing branches, therefore
618 // first we check, if branches with the same title already exist.
619 // If yes - exits without writing.
621 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
623 TClonesArray * trackSegments = gime->TrackSegments() ;
624 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
626 if(fWrite){ //We write TreeT
627 TTree * treeT = gime->TreeT();
630 Int_t bufferSize = 32000 ;
631 TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
634 gime->WriteTracks("OVERWRITE");
635 gime->WriteTrackSegmentMaker("OVERWRITE");
640 //____________________________________________________________________________
641 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
643 // option deb - prints # of found TrackSegments
644 // option deb all - prints as well indexed of found RecParticles assigned to the TS
646 TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ;
648 Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ;
649 printf("nevent: %d\n", gAlice->GetEvNumber()) ;
650 printf(" Found %d TrackSegments\n", trackSegments->GetEntriesFast() );
652 if(strstr(option,"all")) { // printing found TS
653 printf("TrackSegment # EMC RP# CPV RP#\n") ;
655 for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
656 AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ;
657 printf(" %d %d %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ;