]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv1.cxx
Update in the T0 configuration @P2.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSTrackSegmentMakerv1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15 /* $Id$ */
16
17 /* History of cvs commits:
18  *
19  * $Log$
20  * Revision 1.93  2007/10/10 09:05:10  schutz
21  * Changing name QualAss to QA
22  *
23  * Revision 1.92  2007/08/28 12:55:08  policheh
24  * Loaders removed from the reconstruction code (C.Cheshkov)
25  *
26  * Revision 1.91  2007/08/07 14:12:03  kharlov
27  * Quality assurance added (Yves Schutz)
28  *
29  * Revision 1.90  2007/07/11 13:43:30  hristov
30  * New class AliESDEvent, backward compatibility with the old AliESD (Christian)
31  *
32  * Revision 1.89  2007/07/03 08:13:04  kharlov
33  * Bug fix in CPV local coordinates
34  *
35  * Revision 1.88  2007/06/27 09:11:07  kharlov
36  * Bug fix for CPV-EMC distance
37  *
38  * Revision 1.87  2007/05/04 14:49:29  policheh
39  * AliPHOSRecPoint inheritance from AliCluster
40  *
41  * Revision 1.86  2007/04/02 15:00:16  cvetan
42  * No more calls to gAlice in the reconstruction
43  *
44  * Revision 1.85  2007/03/28 19:18:15  kharlov
45  * RecPoints recalculation in TSM removed
46  *
47  * Revision 1.84  2007/03/07 07:01:21  hristov
48  * Fixing copy/paste erro. Additional protections
49  *
50  * Revision 1.83  2007/03/06 21:07:37  kharlov
51  * DP: xz CPV-EMC distance filled to TS
52  *
53  * Revision 1.82  2007/03/06 06:54:48  kharlov
54  * DP:Calculation of cluster properties dep. on vertex added
55  *
56  * Revision 1.81  2007/02/05 10:02:40  kharlov
57  * Module numbering is corrected
58  *
59  * Revision 1.80  2006/08/28 10:01:56  kharlov
60  * Effective C++ warnings fixed (Timur Pocheptsov)
61  *
62  * Revision 1.79  2006/04/25 12:41:15  hristov
63  * Moving non-persistent data to AliESDfriend (Yu.Belikov)
64  *
65  * Revision 1.78  2005/11/18 13:04:51  hristov
66  * Bug fix
67  *
68  * Revision 1.77  2005/11/17 23:34:36  hristov
69  * Corrected logics
70  *
71  * Revision 1.76  2005/11/17 22:29:12  hristov
72  * Faster version, no attempt to match tracks outside the PHOS acceptance
73  *
74  * Revision 1.75  2005/11/17 12:35:27  hristov
75  * Use references instead of objects. Avoid to create objects when they are not really needed
76  *
77  * Revision 1.74  2005/07/08 14:01:36  hristov
78  * Tracking in non-uniform nmagnetic field (Yu.Belikov)
79  *
80  * Revision 1.73  2005/05/28 14:19:05  schutz
81  * Compilation warnings fixed by T.P.
82  *
83  */
84
85 //_________________________________________________________________________
86 // Implementation version 1 of algorithm class to construct PHOS track segments
87 // Track segment for PHOS is list of 
88 //        EMC RecPoint + (possibly) CPV RecPoint
89 // To find TrackSegments we do the following: 
90 //  for each EMC RecPoints we look at
91 //   CPV RecPoints in the radious fRcpv. 
92 //  If there is such a CPV RecPoint, 
93 //   we make "Link" it is just indexes of EMC and CPV RecPoint and distance
94 //   between them in the PHOS plane. 
95 //  Then we sort "Links" and starting from the 
96 //   least "Link" pointing to the unassined EMC and CPV RecPoints assing them to 
97 //   new TrackSegment. 
98 // If there is no CPV RecPoint we make TrackSegment 
99 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
100 //// In principle this class should be called from AliPHOSReconstructor, but 
101 // one can use it as well in standalone mode.
102 // Use  case:
103 //  root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
104 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
105 //               // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
106 //               // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")                       
107 //  root [1] t->ExecuteTask()
108 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
109 //  root [4] t->ExecuteTask("deb all time") 
110 //                 
111 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH) 
112 //
113
114 // --- ROOT system ---
115 #include "AliMagF.h"
116 #include "AliTracker.h"
117 #include "AliGeomManager.h"
118 #include "AliCluster.h"
119 #include "AliKalmanTrack.h"
120 #include "AliGlobalQADataMaker.h"
121
122 #include "TVector3.h"
123 #include "TTree.h"
124 #include "TBenchmark.h"
125
126 // --- Standard library ---
127 #include "Riostream.h"
128 // --- AliRoot header files ---
129 #include "AliPHOSGeometry.h"
130 #include "AliPHOSTrackSegmentMakerv1.h"
131 #include "AliPHOSTrackSegment.h"
132 #include "AliPHOSLink.h"
133 #include "AliESDEvent.h"
134 #include "AliESDtrack.h"
135 #include "AliPHOSEmcRecPoint.h"
136 #include "AliPHOSCpvRecPoint.h"
137 #include "AliLog.h"
138 #include "AliMagF.h"
139
140 ClassImp( AliPHOSTrackSegmentMakerv1) 
141
142
143 //____________________________________________________________________________
144 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() :
145   AliPHOSTrackSegmentMaker(),
146   fDefaultInit(kTRUE),
147   fWrite(kFALSE),
148   fNTrackSegments(0),
149   fRcpv(0.f),
150   fRtpc(0.f),
151   fVtx(0.f), 
152   fLinkUpArray(0),
153   fEmcFirst(0),
154   fEmcLast(0),
155   fCpvFirst(0),
156   fCpvLast(0),
157   fModule(0),
158   fTrackSegments(NULL)
159 {
160   // default ctor (to be used mainly by Streamer)
161   InitParameters() ; 
162 }
163
164 //____________________________________________________________________________
165 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(AliPHOSGeometry *geom) :
166   AliPHOSTrackSegmentMaker(geom),
167   fDefaultInit(kFALSE),
168   fWrite(kFALSE),
169   fNTrackSegments(0),
170   fRcpv(0.f),
171   fRtpc(0.f),
172   fVtx(0.f), 
173   fLinkUpArray(0),
174   fEmcFirst(0),
175   fEmcLast(0),
176   fCpvFirst(0),
177   fCpvLast(0),
178   fModule(0),
179   fTrackSegments(NULL)
180 {
181   // ctor
182   InitParameters() ; 
183   Init() ;
184   fESD = 0;
185 }
186
187
188 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const AliPHOSTrackSegmentMakerv1 & tsm) :
189   AliPHOSTrackSegmentMaker(tsm),
190   fDefaultInit(kFALSE),
191   fWrite(kFALSE),
192   fNTrackSegments(0),
193   fRcpv(0.f),
194   fRtpc(0.f),
195   fVtx(0.f), 
196   fLinkUpArray(0),
197   fEmcFirst(0),
198   fEmcLast(0),
199   fCpvFirst(0),
200   fCpvLast(0),
201   fModule(0),
202   fTrackSegments(NULL)
203 {
204   // cpy ctor: no implementation yet
205   // requested by the Coding Convention
206   Fatal("cpy ctor", "not implemented") ;
207 }
208
209
210 //____________________________________________________________________________
211  AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
212
213   // dtor
214   // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
215   if (!fDefaultInit)  
216     delete fLinkUpArray ;
217   if (fTrackSegments) {
218     fTrackSegments->Delete();
219     delete fTrackSegments;
220   }
221 }
222
223 //____________________________________________________________________________
224 void  AliPHOSTrackSegmentMakerv1::FillOneModule()
225 {
226   // Finds first and last indexes between which 
227   // clusters from one PHOS module are
228
229   //First EMC clusters
230   Int_t totalEmc = fEMCRecPoints->GetEntriesFast() ;
231   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
232         ((dynamic_cast<AliPHOSRecPoint *>(fEMCRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule ); 
233       fEmcLast ++)  ;
234   
235   //Now CPV clusters
236   Int_t totalCpv = fCPVRecPoints->GetEntriesFast() ;
237
238     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) && 
239          ((dynamic_cast<AliPHOSRecPoint *>(fCPVRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule ); 
240        fCpvLast ++) ;
241       
242 }
243
244 //____________________________________________________________________________
245 void  AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
246                                                          AliPHOSCpvRecPoint * cpvClu, 
247                                                          Int_t &trackindex, 
248                                                          Float_t &dx, Float_t &dz) const
249 {
250   // Calculates the distance between the EMC RecPoint and the CPV RecPoint
251   // If no CPV, calculates the distance between the EMC RecPoint and the track
252   // prolongation to the PHOS module plane.
253   // Clusters are sorted in "rows" and "columns" of width 1 cm
254
255 //  Float_t delta = 1 ;  // Width of the rows in sorting of RecPoints (in cm)
256 //                       // if you change this value, change it as well in xxxRecPoint::Compare()
257
258   if(!cpvClu) {
259     
260     trackindex = -1;
261     dx=999.;
262     dz=999.;
263     
264     if(!emcClu) {
265       return;
266     }
267
268     // *** Start the matching
269     Int_t nt=fESD->GetNumberOfTracks();
270     Int_t iPHOSMod = emcClu->GetPHOSMod()  ;
271     //Calculate actual distance to PHOS module
272     TVector3 globaPos ;
273     fGeom->Local2Global(iPHOSMod, 0.,0., globaPos) ;
274     const Double_t rPHOS = globaPos.Pt() ; //Distance to PHOS module
275     const Double_t kYmax = 72.+10. ; //Size of the module (with some reserve) in phi direction
276     const Double_t kZmax = 64.+10. ; //Size of the module (with some reserve) in z direction
277     const Double_t kAlpha= 20./180.*TMath::Pi() ; //TPC sector angular size
278     Double_t minDistance = 1.e6;
279
280     TMatrixF gmat;
281     TVector3 gposRecPoint; // global (in ALICE frame) position of rec. point
282     emcClu->GetGlobalPosition(gposRecPoint,gmat);
283     Double_t gposTrack[3] ; 
284
285     Double_t bz = GetBz() ; //B-Field for approximate matching
286     Double_t b[3]; 
287     for (Int_t i=0; i<nt; i++) {
288       AliESDtrack *esdTrack=fESD->GetTrack(i);
289
290 //     // Skip the tracks having "wrong" status (has to be checked/tuned)
291 //     ULong_t status = esdTrack->GetStatus();
292 //     if ((status & AliESDtrack::kTRDout)   == 0) continue;
293 //     if ((status & AliESDtrack::kTRDrefit) == 1) continue;
294
295       AliExternalTrackParam t(*esdTrack);
296       Int_t isec=Int_t(t.GetAlpha()/kAlpha);
297       Int_t imod=-isec-2; // PHOS module
298
299       Double_t y;                       // Some tracks do not reach the PHOS
300       if (!t.GetYAt(rPHOS,bz,y)) continue; //    because of the bending
301
302       Double_t z; t.GetZAt(rPHOS,bz,z);
303       if (TMath::Abs(z) > kZmax) continue; // Some tracks miss the PHOS in Z
304
305       Bool_t ok=kTRUE;
306       while (TMath::Abs(y) > kYmax) {   // Find the matching module
307         Double_t alp=t.GetAlpha();
308         if (y > kYmax) {
309           if (!t.Rotate(alp+kAlpha)) {ok=kFALSE; break;}
310           imod--;
311         } else if (y < -kYmax) {
312           if (!t.Rotate(alp-kAlpha)) {ok=kFALSE; break;}
313           imod++;
314         }
315         if (!t.GetYAt(rPHOS,bz,y)) {ok=kFALSE; break;}
316       }
317       if (!ok) continue; // Track rotation failed
318
319
320       if(imod!= iPHOSMod-1) 
321         continue; //not even approximate coincidence
322
323       //t.CorrectForMaterial(...); // Correct for the TOF material, if needed
324       t.GetBxByBz(b) ;
325       t.PropagateToBxByBz(rPHOS,b);        // Propagate to the matching module
326       t.GetXYZ(gposTrack) ;
327       
328       Double_t ddx = gposTrack[0] - gposRecPoint.X(), ddy = gposTrack[1] - gposRecPoint.Y(), ddz = gposTrack[2] - gposRecPoint.Z();
329       Double_t d2 = ddx*ddx + ddy*ddy + ddz*ddz;
330       if(d2 < minDistance) {
331         dx = TMath::Sign(TMath::Sqrt(ddx*ddx + ddy*ddy),ddx) ;
332         dz = ddz ;
333         trackindex=i;
334         minDistance=d2 ;
335       }
336     }
337     return ;
338   }
339
340
341   TVector3 emcGlobal;
342   fGeom->GetGlobalPHOS((AliPHOSRecPoint*)emcClu,emcGlobal);
343     
344   // Radius from IP to current point
345   Double_t rEMC = TMath::Abs(emcGlobal.Pt());
346     
347   // Extrapolate the global track direction to EMC
348   // and find the closest track
349     
350   Int_t nTracks = fESD->GetNumberOfTracks();
351     
352   AliESDtrack *track;
353   Double_t xyz[] = {-1,-1,-1};
354   Double_t pxyz[3];
355   Double_t zEMC,xEMC;
356   Int_t module;
357   TVector3 vecP;
358   TVector3 locClu;
359
360   Float_t minDistance = 1.e6;
361   Float_t dr;
362
363   for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
364     track = fESD->GetTrack(iTrack);
365     if (!track->GetXYZAt(rEMC, fESD->GetMagneticField(), xyz)) continue;
366     
367     AliDebug(1,Form("Event %d, iTrack: %d, (%.3f,%.3f,%.3f)",
368                     fESD->GetEventNumberInFile(),iTrack,xyz[0],xyz[1],xyz[2]));
369       
370     if (track->GetPxPyPzAt(rEMC,fESD->GetMagneticField(),pxyz)) { 
371
372       vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
373       fGeom->ImpactOnEmc(xyz,vecP.Theta(),vecP.Phi(),module,zEMC,xEMC) ;
374
375       if(!module) continue;
376       AliDebug(1,Form("\t\tTrack hit PHOS! Module: %d, (x,z)=(%.3f,%.3f)",module,xEMC,zEMC));
377         
378       if(emcClu->GetPHOSMod() != module) continue;
379         
380       // match track to EMC cluster
381       emcClu->GetLocalPosition(locClu);
382         
383       Float_t delta_x = xEMC - locClu.X();
384       Float_t delta_z = zEMC - locClu.Z();
385       dr = TMath::Sqrt(delta_x*delta_x + delta_z*delta_z);
386       AliDebug(1,Form("\tMatch iTrack=%d: (dx,dz)=(%.3f,%.3f)",iTrack,delta_x,delta_z));
387         
388       if(dr<minDistance) {
389         trackindex = iTrack;
390         minDistance = dr;
391         dx = delta_x;
392         dz = delta_z;
393       }
394     }
395       
396   }
397     
398   if(trackindex>=0)
399     AliDebug(1,Form("\t\tBest match for (xClu,zClu,eClu)=(%.3f,%.3f,%.3f): iTrack=%d, dR=%.3f",
400                     locClu.X(),locClu.Z(),emcClu->GetEnergy(),
401                     trackindex,TMath::Sqrt(dx*dx+dz*dz)));        
402   return;
403   
404   Float_t distance2Track = fRtpc ; 
405   
406   trackindex = -1 ; // closest track within fRCpv 
407   
408   TVector3 vecEmc ;   // Local position of EMC recpoint
409   TVector3 vecPloc ;     // Momentum direction at CPV plain
410   
411   //toofar = kTRUE ;
412   if(emcClu->GetPHOSMod() != cpvClu->GetPHOSMod()){
413     dx=999. ;
414     dz=999. ;
415     return ;
416   }
417   
418   emcClu->GetLocalPosition(vecEmc) ;
419   
420   Double_t xCPV,zCPV ; //EMC-projected coordinates of CPV cluster 
421   TVector3 cpvGlobal; // Global position of the CPV recpoint
422   fGeom->GetGlobalPHOS((AliPHOSRecPoint*)cpvClu,cpvGlobal);
423   Double_t vtxCPV[3]={cpvGlobal.X(),cpvGlobal.Y(),cpvGlobal.Z()} ;
424   Int_t dummyMod ;
425
426   if (fESD == 0x0) {
427     //if no track information available, assume straight line from IP to emcal
428     fGeom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,zCPV,xCPV) ;
429     dx=xCPV - vecEmc.X() ;
430     dz=zCPV - vecEmc.Z() ;
431     return ;
432   } 
433   
434   //if there is ESD try to correct distance using TPC information on particle direct in CPV 
435   if (fESD != 0x0) {
436
437     Double_t rCPV = cpvGlobal.Pt() ;// Radius from IP to current point 
438
439     // Extrapolate the global track direction if any to CPV and find the closest track
440     Int_t nTracks = fESD->GetNumberOfTracks();
441     Int_t iClosestTrack = -1;
442     Double_t minDistance = 1.e6;
443     TVector3 inPHOS ; 
444
445     AliESDtrack *track;
446     Double_t xyz[3] ;
447     Double_t pxyz[3]; 
448     for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
449       track = fESD->GetTrack(iTrack);
450       if (!track->GetXYZAt(rCPV, fESD->GetMagneticField(), xyz))
451         continue; //track coord on the cylinder of PHOS radius
452       if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
453         continue;
454       //Check if this track hits PHOS
455       inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
456       distance2Track = inPHOS.Angle(cpvGlobal) ;
457       // Find the closest track to the CPV recpoint
458       if (distance2Track < minDistance) {
459         minDistance = distance2Track;
460         iClosestTrack = iTrack;
461       }
462     }
463
464     if (iClosestTrack != -1) {
465       track = fESD->GetTrack(iClosestTrack);
466       if (track->GetPxPyPzAt(rCPV, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
467         vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
468         fGeom->ImpactOnEmc(vtxCPV,vecP.Theta(),vecP.Phi(),dummyMod,zCPV,xCPV) ;
469       }
470     }
471     
472     if(minDistance < fRtpc ){
473       trackindex = iClosestTrack ; 
474     }
475   }
476   if(trackindex!=-1){
477     // If the closest global track is found, calculate EMC-CPV distance from it
478     dx=xCPV - vecEmc.X() ;
479     dz=zCPV - vecEmc.Z() ;
480   }
481   else{
482     // If no global track was found, just take the nearest CPV point
483     fGeom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,zCPV,xCPV) ;
484     dx=xCPV - vecEmc.X() ;
485     dz=zCPV - vecEmc.Z() ;
486   }
487   return ;
488 }
489 //____________________________________________________________________________
490 void  AliPHOSTrackSegmentMakerv1::Init()
491 {
492   // Make all memory allocations that are not possible in default constructor
493   
494   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
495   fTrackSegments = new TClonesArray("AliPHOSTrackSegment",100);
496   fTrackSegments->SetName("TRACKS");
497 }
498
499 //____________________________________________________________________________
500 void  AliPHOSTrackSegmentMakerv1::InitParameters()
501 {
502   //Initializes parameters
503   fRcpv      = 10. ;
504   fRtpc      = 4. ;
505   fEmcFirst  = 0 ;    
506   fEmcLast   = 0 ;   
507   fCpvFirst  = 0 ;   
508   fCpvLast   = 0 ;   
509   fLinkUpArray = 0 ;
510   fWrite                   = kTRUE ;
511 }
512
513
514 //____________________________________________________________________________
515 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
516
517   // Finds distances (links) between all EMC and CPV clusters, 
518   // which are not further apart from each other than fRcpv 
519   // and sort them in accordance with this distance
520   
521   fLinkUpArray->Clear() ;    
522
523   AliPHOSCpvRecPoint * cpv ;
524   AliPHOSEmcRecPoint * emcclu ;
525
526   Int_t iLinkUp  = 0 ;
527   
528   Int_t iEmcRP;
529   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
530     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(fEMCRecPoints->At(iEmcRP)) ;
531
532     //Bool_t toofar ;        
533     Int_t iCpv = 0 ;    
534     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
535       
536       cpv = dynamic_cast<AliPHOSCpvRecPoint *>(fCPVRecPoints->At(iCpv)) ;
537       Int_t track = -1 ; 
538       Float_t dx,dz ;
539       GetDistanceInPHOSPlane(emcclu, cpv, track,dx,dz) ;     
540       if(TMath::Sqrt(dx*dx+dz*dz) < fRcpv ){ 
541         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(dx, dz, iEmcRP, iCpv, track) ;
542       }      
543     }
544   } 
545   
546   fLinkUpArray->Sort() ;  //first links with smallest distances
547 }
548
549 //____________________________________________________________________________
550 void  AliPHOSTrackSegmentMakerv1::MakePairs()
551
552   // Using the previously made list of "links", we found the smallest link - i.e. 
553   // link with the least distance between EMC and CPV and pointing to still 
554   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
555   // remove them from the list of "unassigned". 
556
557   //Make arrays to mark clusters already chosen
558   Int_t * emcExist = 0;
559   if(fEmcLast > fEmcFirst)
560     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
561   
562   Int_t index;
563   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
564     emcExist[index] = 1 ;
565   
566   Bool_t * cpvExist = 0;
567   if(fCpvLast > fCpvFirst)
568     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
569   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
570     cpvExist[index] = kTRUE ;
571   
572   
573   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
574   TIter nextUp(fLinkUpArray) ;
575   
576   AliPHOSLink * linkUp ;
577   
578   AliPHOSCpvRecPoint * nullpointer = 0 ;
579   
580   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
581
582     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
583
584       if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
585          Float_t dx,dz ;
586          linkUp->GetXZ(dx,dz) ;
587          new ((* fTrackSegments)[fNTrackSegments]) 
588            AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(fEMCRecPoints->At(linkUp->GetEmc())) , 
589                                dynamic_cast<AliPHOSCpvRecPoint *>(fCPVRecPoints->At(linkUp->GetCpv())) , 
590                                linkUp->GetTrack(),dx,dz) ;
591          
592        (dynamic_cast<AliPHOSTrackSegment *>(fTrackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
593        fNTrackSegments++ ;
594        emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
595        //mark CPV recpoint as already used 
596        cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
597       } //if CpvUp still exist
598     } 
599   }        
600
601   //look through emc recPoints left without CPV
602   if(emcExist){ //if there is emc rec point
603     Int_t iEmcRP ;
604     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
605       if(emcExist[iEmcRP] > 0 ){
606         Int_t track = -1 ;
607         Float_t dx,dz ;
608         AliPHOSEmcRecPoint *emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(fEMCRecPoints->At(iEmcRP+fEmcFirst));
609         GetDistanceInPHOSPlane(emcclu, 0, track,dx,dz);
610         if(track<0)
611           new ((*fTrackSegments)[fNTrackSegments]) AliPHOSTrackSegment(emcclu,nullpointer) ;
612         else
613           new ((*fTrackSegments)[fNTrackSegments]) AliPHOSTrackSegment(emcclu,0,track,dx,dz);
614         (dynamic_cast<AliPHOSTrackSegment *>(fTrackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
615         fNTrackSegments++;    
616       } 
617     }
618   }
619   delete [] emcExist ; 
620   delete [] cpvExist ; 
621 }
622
623 //____________________________________________________________________________
624 void AliPHOSTrackSegmentMakerv1::Clusters2TrackSegments(Option_t *option)
625 {
626   // Steering method to perform track segment construction for the current event
627   // Returns an array with the found track-segments.
628   
629   if(strstr(option,"tim"))
630     gBenchmark->Start("PHOSTSMaker");
631  
632   if(strstr(option,"print")) {
633     Print() ; 
634     return ; 
635   }
636   
637   //Make some initializations 
638   fNTrackSegments = 0 ;
639   fEmcFirst = 0 ;    
640   fEmcLast  = 0 ;   
641   fCpvFirst = 0 ;   
642   fCpvLast  = 0 ;   
643
644   fTrackSegments->Clear();
645
646   //   if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
647
648   for(fModule = 1; fModule <= fGeom->GetNModules() ; fModule++ ) {
649     FillOneModule() ; 
650     MakeLinks() ;
651     MakePairs() ;
652   }
653     
654   if(strstr(option,"deb"))
655     PrintTrackSegments(option);
656
657   if(strstr(option,"tim")){
658     gBenchmark->Stop("PHOSTSMaker");
659     Info("Exec", "took %f seconds for making TS", 
660          gBenchmark->GetCpuTime("PHOSTSMaker")); 
661   }
662 }
663
664 //____________________________________________________________________________
665 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
666 {
667   //  Print TrackSegmentMaker parameters
668
669   TString message("") ;
670   if( strcmp(GetName(), "") != 0 ) {
671     message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ; 
672     message += "Making Track segments\n" ;
673     message += "with parameters:\n" ; 
674     message += "     Maximal EMC - CPV distance (cm) %f\n" ;
675     message += "============================================\n" ;
676     Info("Print", message.Data(),fRcpv) ;
677   }
678   else
679     Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
680 }
681
682 //____________________________________________________________________________
683 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
684 {
685   // option deb - prints # of found TrackSegments
686   // option deb all - prints as well indexed of found RecParticles assigned to the TS
687
688   Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ; 
689   printf("        Found %d TrackSegments\n", fTrackSegments->GetEntriesFast() ); 
690   
691   if(strstr(option,"all")) {  // printing found TS
692     printf("TrackSegment #  EMC RP#  CPV RP#\n") ; 
693     Int_t index;
694     for (index = 0 ; index <fTrackSegments->GetEntriesFast() ; index++) {
695       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )fTrackSegments->At(index) ; 
696       printf("   %d           %d        %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ; 
697     }   
698   }
699 }
700 //__________________________________________________________________________
701 Double_t AliPHOSTrackSegmentMakerv1::GetBz()const
702 {
703   Double_t kAlmost0Field=1.e-13;
704   AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
705   if (!fld) return 0.5*kAlmost0Field;
706   Double_t bz = fld->SolenoidField();
707   return TMath::Sign(0.5*kAlmost0Field,bz) + bz;
708 }
709 //__________________________________________________________________________
710 void AliPHOSTrackSegmentMakerv1::GetBxByBz(const Double_t r[3], Double_t b[3])const {
711   //------------------------------------------------------------------
712   // Returns Bx, By and Bz (kG) at the point "r" .
713   //------------------------------------------------------------------
714   Double_t kAlmost0Field=1.e-13;
715   AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
716   if (!fld) {
717      b[0] = b[1] = 0.;
718      b[2] = 0.5*kAlmost0Field;
719      return;
720   }
721
722   if (fld->IsUniform()) {
723      b[0] = b[1] = 0.;
724      b[2] = fld->SolenoidField();
725   }  else {
726      fld->Field(r,b);
727   }
728   b[2] = (TMath::Sign(0.5*kAlmost0Field,b[2]) + b[2]);
729   return;
730 }
731
732