]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDTracker.cxx
Fixed event selection bug for pPb
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDTracker.cxx
1 #include "AliHMPIDTracker.h"     //class header
2 #include "AliHMPIDtrack.h"       //class header
3 #include "AliHMPIDCluster.h"     //GetTrackPoint(),PropagateBack() 
4 #include "AliHMPIDParam.h"       //GetTrackPoint(),PropagateBack()
5 #include "AliHMPIDPid.h"         //Recon(),reconHTA()
6 #include "AliHMPIDRecon.h"       //Recon()
7 #include "AliHMPIDRecoParamV1.h"   //Recon()
8 #include "AliHMPIDReconstructor.h"//Recon()
9 #include "AliHMPIDReconHTA.h"    //ReconHTA()
10 #include <AliLog.h>              //Recon()  
11 #include <AliESDEvent.h>         //PropagateBack(),Recon()  
12 #include <AliESDtrack.h>         //Intersect() 
13 #include <AliTracker.h> 
14 #include <AliRun.h>              //GetTrackPoint(),PropagateBack()  
15 #include <AliTrackPointArray.h>  //GetTrackPoint()
16 #include <AliAlignObj.h>         //GetTrackPoint()
17 #include <AliCDBManager.h>       //PropageteBack()
18 #include <AliCDBEntry.h>         //PropageteBack()
19 #include "TTreeStream.h"         // debug streamer
20 //
21 // HMPID base class fo tracking
22 //.
23 //.
24 //.
25 ClassImp(AliHMPIDTracker)
26 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 AliHMPIDTracker::AliHMPIDTracker():
28   AliTracker(),
29   fClu(new TObjArray(AliHMPIDParam::kMaxCh+1)),
30   fDebugStreamer(0)
31 {
32 // ctor. Create TObjArray of TClonesArray of AliHMPIDCluster  
33 // 
34 //  
35   fClu->SetOwner(kTRUE);
36   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
37   fDebugStreamer = new TTreeSRedirector("HMPIDdebug.root");
38
39 }//ctor
40
41
42 AliHMPIDTracker::~AliHMPIDTracker(){
43   //
44   // destructor
45   // 
46   delete fClu;
47   if (fDebugStreamer) delete fDebugStreamer;
48 }
49
50 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++    
51 Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
52 {
53 // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
54 // MIP cluster is reffered by index which is stored in AliESDtrack  ???????
55 // Arguments: idx- cluster index which is stored by HMPID in AliESDtrack
56 //            point- reference to the object where to store the point     
57 //   Returns: status of operation  if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. 
58   if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
59   Int_t iCham=idx/1000000; Int_t iClu=idx%1000000;
60   iClu = iClu%1000; //GetHMPIDcluIdx -> 1e+6*ch + 1e+3*clusize + cluIdx;
61   point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham));//layer and chamber number
62   TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham];
63   if(!pArr) return kFALSE;
64   AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
65   if(!pClu) return kFALSE;
66   Float_t xyz[3];
67   pClu->GetGlobalXYZ(xyz);
68   Float_t cov[6];
69   pClu->GetGlobalCov(cov);
70   point.SetXYZ(xyz,cov);
71   return kTRUE;
72 }//GetTrackPoint()
73 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
74 Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
75 {
76 // Static method to find intersection in between given track and HMPID chambers
77 // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
78 //   Returns: intersected chamber ID or -1
79   AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
80   for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){                              //chambers loop
81     Int_t chInt = IntTrkCha(i,hmpTrk,xPc,yPc,xRa,yRa,theta,phi);
82     if(chInt>=0) {delete hmpTrk;hmpTrk=0x0;return chInt;}
83   }                                                                                             //chambers loop
84   delete hmpTrk; hmpTrk=0x0;
85   return -1;                                                                                    //no intersection with HMPID chambers
86 }//IntTrkCha()
87 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
88 Int_t AliHMPIDTracker::IntTrkCha(Int_t ch,AliHMPIDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
89 {
90 // Static method to find intersection in between given track and HMPID chambers
91 // Arguments: pTrk- HMPID track; xPc,yPc- track intersection with PC in LORS [cm]
92 //   Returns: intersected chamber ID or -1
93     AliHMPIDParam *pParam=AliHMPIDParam::Instance();
94     Double_t p1[3],n1[3]; 
95     pParam->Norm(ch,n1); 
96     pParam->Point(ch,p1,AliHMPIDParam::kRad);                                                    //point & norm  for middle of radiator plane
97     Double_t p2[3],n2[3]; 
98     pParam->Norm(ch,n2);
99     pParam->Point(ch,p2,AliHMPIDParam::kPc);                                                     //point & norm  for entrance to PC plane
100     if(pTrk->Intersect(p1,n1)==kFALSE) return -1;                                                //try to intersect track with the middle of radiator
101     if(pTrk->Intersect(p2,n2)==kFALSE) return -1;   
102     pParam->Mars2LorsVec(ch,n1,theta,phi);                                                       //track angles at RAD
103     pParam->Mars2Lors   (ch,p1,xRa,yRa);                                                         //TRKxRAD position
104     pParam->Mars2Lors   (ch,p2,xPc,yPc);                                                         //TRKxPC position
105     if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return ch;                     //return intersected chamber  
106   return -1;                                                                                     //no intersection with HMPID chambers
107 }//IntTrkCha()
108 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
109 Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
110 {
111 // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
112 // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
113 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
114   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
115   pCluTree->GetEntry(0);
116   return 0;  
117 }//LoadClusters()
118 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
119 Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
120 {
121 // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
122 // Agruments: pEsd - pointer to ESD
123 //   Returns: error code    
124   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
125   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
126   if(!pNmeanEnt) AliError("No Nmean C6F14 ");
127   if(!pQthreEnt) AliError("No Qthre");
128     
129   return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());  
130 }//PropagateBack()
131 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
132 Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
133 {
134 // Static method to reconstruct the Cherenkov angle for all valid tracks of a given event.
135 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
136 // Returns: error code, 0 if no errors   
137 //
138 // Algortihm: Loop over tracks
139 //
140 // 1. Find the closest MIP cluster using fast tracks extrapolation method
141 // 2. Propagate track to the MIP cluster using the STEER method
142 // 3. Update the track information with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
143 // 4. Propagate back the constrained track to the radiator radius ( not exact yet)
144 // 5. Propagation in the last 10 cm with the fast method 
145 // 6. Set ESDtrack information
146 // 7. Calculate the Cherenkov angle
147
148   const Double_t kMaxSnp=0.9;   //maximal snp for prolongation
149   const Double_t kdRadiator=10; // distance between radiator and the plane
150   AliHMPIDRecon recon;                                                                           //instance of reconstruction class, nothing important in ctor
151   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
152   Float_t xPc,yPc,xRa,yRa,theta,phi;
153   
154   Double_t cluLORS[2]={0};
155   Int_t nMipClusTot=0;
156   
157   Double_t qthre = 0;   Double_t nmean=0; Int_t hvsec=0;
158   Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
159
160   Bool_t tsRight = kTRUE;
161   
162   UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin();                                        //
163   UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax();                                        //
164   UInt_t ts = pEsd->GetTimeStamp();
165   
166   if(ts<tsmin || ts>tsmax) {
167     AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
168     tsRight = kFALSE;
169   }
170    
171   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
172            
173     //    Double_t bestChi2=99999;chi2=99999;                                                     //init. track matching params
174     Double_t dmin=999999,bz=0,distCut=1,distParams[5]={1};
175
176     Bool_t isOkDcut=kFALSE;
177     Bool_t isOkQcut=kFALSE;
178     Bool_t isMatched=kFALSE;
179     
180     AliHMPIDCluster *bestHmpCluster=0x0;                                                          //the best matching cluster
181     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //get reconstructed track   
182     
183     if(!pTrk->IsOn(AliESDtrack::kTPCout)) continue;
184  
185     if(pTrk->IsOn(AliESDtrack::kTPCrefit)) continue;
186     AliESDfriendTrack *ftrack= (AliESDfriendTrack *)pTrk->GetFriendTrack();
187     if (!ftrack) continue; 
188     if (!ftrack->GetTPCOut()) continue;
189     AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
190     AliHMPIDtrack *hmpTrkConstrained = 0;                                                         //create a hmpid track to be used for propagation and matching 
191     hmpTrk->Set(ftrack->GetTPCOut()->GetX(), ftrack->GetTPCOut()->GetAlpha(),ftrack->GetTPCOut()->GetParameter(), ftrack->GetTPCOut()->GetCovariance());
192     //
193     bz=AliTracker::GetBz();
194           
195     //initial flags for HMPID ESD infos    
196     pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
197     pTrk->SetHMPIDmip(0,0,0,0);                                                                //store mip info in any case 
198     pTrk->SetHMPIDcluIdx(99,99999);                                                            //chamber not found, mip not yet considered
199     pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
200     
201     Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);                                        //find the intersected chamber for this track 
202     if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;}                                                                         //no intersection at all, go after next track
203
204     pTrk->SetHMPIDtrk(xPc,yPc,theta,phi);                                                        //store initial infos
205     pTrk->SetHMPIDcluIdx(ipCh,9999);                                                             //set chamber, index of cluster + cluster size
206     
207     // track intersects the chamber ipCh: find the MIP          
208     
209     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
210     nMipClusTot = pMipCluLst->GetEntries();                                                     //total number of clusters in the given chamber
211     nClusCh[ipCh] = nMipClusTot;
212     
213     if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;}                                                                         
214     
215     Int_t index=-1;                                                                               //index of the "best" matching cluster
216     //
217     // 1. Find the closest MIP cluster using fast tracks extrapolation method
218     //
219     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                                 //clusters loop
220       
221       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                      //get the cluster
222       // evaluate qThre
223       if(tsRight){
224         if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {
225           qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts);                                          //
226         } else {                                                                                   // in the past just 1 qthre
227           hvsec = pParam->InHVSector(pClu->Y());                                                   //  per chamber
228           if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts);                                     //
229         } 
230       } else qthre = pParam->QCut();
231       
232       //      
233       if(pClu->Q()<qthre) continue;                                                                      //charge compartible with MIP clusters      
234       isOkQcut = kTRUE;
235       //
236       cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y();                                                        //get the LORS coordinates of the cluster
237       Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
238       
239       if(dist<dmin) {
240         dmin = dist;
241         index=iClu;
242         bestHmpCluster=pClu;
243       }
244     } // clusters loop
245
246     // moved down
247     /*if(!isOkQcut) {
248       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
249       delete hmpTrk;hmpTrk=0x0; continue;                                                                     
250     }*/
251     
252     //
253     // 2. Propagate track to the MIP cluster using the STEER method
254     //
255         
256     if(!bestHmpCluster) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0;  continue;}     
257         
258     TVector3 vG = pParam->Lors2Mars(ipCh,bestHmpCluster->X(),bestHmpCluster->Y());    
259     Double_t gx = vG[0];
260     Double_t gy = vG[1];
261     Double_t gz = vG[2];
262     Double_t alpha=TMath::ATan2(gy,gx);
263     Double_t radiusH=TMath::Sqrt(gy*gy+gx*gx);
264     if (!(hmpTrk->Rotate(alpha,kTRUE))) continue;
265     if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrk,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}    
266     //
267     // 3. Update the track with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
268     //
269     AliExternalTrackParam trackC(*hmpTrk);
270     Double_t posC[2]={0,gz};
271     Double_t covC[3]={0.1*0.1, 0, 0.1*0.1};
272     trackC.Update(posC,covC);
273     //
274     // 4. Propagate back the constrained track to the radiator radius ( not exact yet)
275     //
276     hmpTrkConstrained = new AliHMPIDtrack(*pTrk);
277     hmpTrkConstrained->Set(trackC.GetX(), trackC.GetAlpha(),trackC.GetParameter(), trackC.GetCovariance());
278     if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrkConstrained,radiusH-kdRadiator,pTrk->GetMass(),1,kFALSE,kMaxSnp,1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;}
279     //
280     // 5. Propagation in the last 10 cm with the fast method 
281     //
282     Float_t xPc0=0, yPc0=0;
283     
284     IntTrkCha(ipCh, hmpTrk, xPc0,yPc0,xRa,yRa,theta,phi);  
285     IntTrkCha(ipCh, hmpTrkConstrained, xPc,yPc,xRa,yRa,theta,phi);  
286     //
287     // 6. Set ESDtrack information
288     //
289     Int_t cluSiz = bestHmpCluster->Size();
290     pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0);  //store mip info in any case 
291     pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                             //set chamber, index of cluster + cluster size    
292     pTrk->SetHMPIDtrk(xPc0,yPc0,theta,phi);
293     //
294     //
295     // Dump debug info if specified
296     // 
297     if (AliHMPIDReconstructor::StreamLevel()>0) {       
298       AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
299       AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);      
300       if(!trackTPC->Rotate(alpha)) continue;
301       if(!trackCurrent->Rotate(alpha)) continue;      
302       Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);      
303       Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);  
304       AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
305       if(!trackTPCNB->Rotate(alpha)) continue;
306       Bool_t statusTPCNB=kTRUE;
307       Double_t bfield[3]={0,0,0};
308       for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){
309         Double_t xyz[3];
310         trackTPCNB->GetXYZ(xyz);
311         GetBxByBz(xyz,bfield);
312         statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield);
313       }
314       statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield);
315       //        
316       Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
317       Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());    
318       //
319       AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
320       Double_t pos[2]={0,gz};
321       Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
322       Double_t chi2C =  trackTPCConstrained->GetPredictedChi2(pos,cov);
323       trackTPCConstrained->Update(pos,cov);
324       (*fDebugStreamer)<<"track"<<
325         "rH="<<radiusH<<                      // radius of cluster
326         "angle="<<tanAlpha<<                  // tan of the local inlination angle
327         "dC="<<deltaC<<                       // delta of the curvature
328         "trackTPC.="<<trackTPC<<              // TPC outer param extrapolated to the HMPID
329         "trackTPCNB.="<<trackTPCNB<<          // TPC track prpagated with material budget correction
330         "chi2C="<<chi2C<<
331         "trackTPCC.="<<trackTPCConstrained<<  // TPC outer param extrapolated to the HMPID constrained
332         "trackCurrent.="<<trackCurrent<<      // current track extrapolated to the HMPID
333         "sTPC="<<statusTPC<<                  // status for the current TPC  track
334         "sCurrent="<<statusCurrent<<          // status for the current global track
335         "cl.="<<bestHmpCluster<<              // HMPID cluster
336         //
337         "t.="<<pTrk<<                        // curent esd track
338         "ft.="<<ftrack<<                     // friend track
339         "hmpTrk.="<<hmpTrk<<                 // hmpid tracks as used in the following code
340         "hmpTrkC.="<<hmpTrkConstrained<<     // constrained hmpid tracks as used in the following code
341         "gx="<<gx<<                          // global cluster position X
342         "gy="<<gy<<                          // Y
343         "gz="<<gz<<                          // Z
344         "\n";
345     }                 
346     //
347     //
348     //    
349     if(!isOkQcut) {
350       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
351       delete hmpTrk;hmpTrk=0x0; 
352       delete hmpTrkConstrained;hmpTrkConstrained=0x0; 
353       continue;                                                                     
354     }        
355     
356     if(AliHMPIDReconstructor::GetRecoParam())                                                 //retrieve distance cut
357     {
358       if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE)                      //distance cut is fixed number
359       { 
360         distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist();
361       }
362       else 
363       {
364         for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar);      //prevision: distance cut is function of momentum
365         distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise
366       }
367     }
368     else {distCut=pParam->DistCut();}
369      
370     //dmin recalculated
371     
372     dmin = TMath::Sqrt((xPc0-bestHmpCluster->X())*(xPc0-bestHmpCluster->X())+(yPc0-bestHmpCluster->Y())*(yPc0-bestHmpCluster->Y()));
373              
374     if(dmin < distCut) {
375       isOkDcut = kTRUE;
376     }   
377     //isOkDcut = kTRUE; // switch OFF cut
378     
379     if(!isOkDcut) {
380       pTrk->SetHMPIDsignal(pParam->kMipDistCut);                                                //closest cluster with enough charge is still too far from intersection
381     }
382     
383     if(isOkQcut*isOkDcut) isMatched = kTRUE;                                                    // MIP-Track matched !!    
384     
385     if(!isMatched) {delete hmpTrk;hmpTrk=0x0;delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;}                                           // If matched continue...
386
387     Bool_t isOk = kTRUE; 
388     if(!isOk) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}
389     pTrk->SetOuterHmpParam(hmpTrkConstrained,AliESDtrack::kHMPIDout);                 
390
391     FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
392  
393     Int_t iRad = pParam->Radiator(yRa);                                                          //evaluate the radiator involved
394     
395     //evaluate nMean
396     if(tsRight){
397      if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
398        nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts);                                               //C6F14 Nmean for this chamber
399      } else {
400        if(iRad < 0) {
401         nmean = -1;
402        } else {
403        Double_t tLow  = ((TF1*)pNmean->At(6*ipCh+2*iRad  ))->Eval(ts);                           //C6F14 low  temp for this chamber
404        Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts);                           //C6F14 high temp for this chamber
405        Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
406        nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
407        }
408        if(nmean < 0){                                                                            //track didn' t pass through the radiator
409          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
410          pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                          //set index of cluster
411          delete hmpTrk;hmpTrk=0x0; 
412          delete hmpTrkConstrained;hmpTrkConstrained=0x0; 
413          continue;
414         }
415      }
416     } else nmean = pParam->MeanIdxRad();    
417     //
418     // 7. Calculate the Cherenkov angle
419     //
420     recon.SetImpPC(xPc0,yPc0);                                                                     //store track impact to PC           
421     recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa);                   //search for Cerenkov angle of this track
422     
423     Double_t thetaCkov = pTrk->GetHMPIDsignal();
424     
425     if (AliHMPIDReconstructor::StreamLevel()>0) {       
426       AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
427       AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);      
428       if(!trackTPC->Rotate(alpha)) continue;
429       if(!trackCurrent->Rotate(alpha)) continue;      
430       Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);      
431       Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);    
432       Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
433       Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());    
434       //
435       AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
436       if(!trackTPCNB->Rotate(alpha)) continue;
437       Bool_t statusTPCNB=kTRUE;
438       Double_t bfield[3]={0,0,0};
439       for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){
440         Double_t xyz[3];
441         trackTPCNB->GetXYZ(xyz);
442         GetBxByBz(xyz,bfield);
443         statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield);
444       }
445       statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield);
446
447       AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
448       Double_t pos[2]={0,gz};
449       Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
450       Double_t chi2C =  trackTPCConstrained->GetPredictedChi2(pos,cov);
451       trackTPCConstrained->Update(pos,cov);
452       (*fDebugStreamer)<<"track2"<<
453         "rH="<<radiusH<<                      // radius of cluster
454         "angle="<<tanAlpha<<                  // tan of the local inlination angle
455         "dC="<<deltaC<<                       // delta of the curvature
456         "trackTPC.="<<trackTPC<<              // TPC outer param extrapolated to the HMPID
457         "trackTPCNB.="<<trackTPCNB<<          // TPC outer param extrapolated to the HMPID
458         "chi2C="<<chi2C<<
459         "trackTPCC.="<<trackTPCConstrained<<  // TPC outer param extrapolated to the HMPID constrained
460         "trackCurrent.="<<trackCurrent<<      // current track extrapolated to the HMPID
461         "sTPC="<<statusTPC<<                  // status for the current TPC  track
462         "sCurrent="<<statusCurrent<<          // status for the current global track
463         "cl.="<<bestHmpCluster<<              // HMPID cluster
464         //
465         "t.="<<pTrk<<                         // curent esd track
466         "ft.="<<ftrack<<                      // friend track
467         "hmpTrk.="<<hmpTrk<<                  // hmpid tracks as used in the following code
468         "hmpTrkC.="<<hmpTrkConstrained<<      // constrained hmpid tracks as used in the following code
469         "gx="<<gx<<                           // global cluster position X
470         "gy="<<gy<<                           // Y
471         "gz="<<gz<<                           // Z
472         "thetaCkov="<<thetaCkov<<             // Cherenkov angle
473         "\n";
474     } 
475                 
476     if(pTrk->GetHMPIDsignal()<0) {
477       delete hmpTrk;hmpTrk=0x0;
478       delete hmpTrkConstrained; hmpTrkConstrained=0x0;
479       continue;}
480         
481     AliHMPIDPid pID;
482     Double_t prob[5];
483     pID.FindPid(pTrk,nmean,5,prob);
484     pTrk->SetHMPIDpid(prob);
485     delete hmpTrk; hmpTrk=0x0;
486     delete hmpTrkConstrained; hmpTrkConstrained=0x0;
487   }//iTrk
488
489   return 0; // error code: 0=no error;
490 }//Recon()
491 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
492 Int_t AliHMPIDTracker::ReconFromKin(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
493 {
494 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
495 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
496 //   Returns: error code, 0 if no errors   
497   
498   AliHMPIDRecon recon;                                                                           //instance of reconstruction class, nothing important in ctor
499   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
500   Float_t xPc,yPc,xRa,yRa,theta,phi;
501   Double_t cluLORS[2]={0};
502   Int_t nMipClusTot=0;
503   Double_t qthre = 0;   Double_t nmean=0; Int_t hvsec=0;
504   Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
505
506   Bool_t tsRight = kTRUE;
507   
508   UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin();                                        //
509   UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax();                                        //
510   UInt_t ts = pEsd->GetTimeStamp();
511   
512   if(ts<tsmin || ts>tsmax) {
513     AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
514     tsRight = kFALSE;
515   }
516    
517   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
518            
519     Double_t dmin=999999,bz=0,distCut=1,distParams[5]={1};
520
521     Bool_t isOkDcut=kFALSE;
522     Bool_t isOkQcut=kFALSE;
523     Bool_t isMatched=kFALSE;
524     
525     AliHMPIDCluster *bestHmpCluster=0x0;                                                          //the best matching cluster
526     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //get reconstructed track   
527     
528     AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
529     bz=AliTracker::GetBz();  
530 //initial flags for HMPID ESD infos    
531     pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
532     pTrk->SetHMPIDmip(0,0,0,0);                                                                //store mip info in any case 
533     pTrk->SetHMPIDcluIdx(99,99999);                                                            //chamber not found, mip not yet considered
534     pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
535     
536     Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);                                        //find the intersected chamber for this track 
537     if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;}                                                                         //no intersection at all, go after next track
538
539     pTrk->SetHMPIDtrk(xPc,yPc,theta,phi);                                                        //store initial infos
540     pTrk->SetHMPIDcluIdx(ipCh,9999);                                                             //set chamber, index of cluster + cluster size
541     
542 // track intersects the chamber ipCh: find the MIP          
543     
544     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
545     nMipClusTot = pMipCluLst->GetEntries();                                                     //total number of clusters in the given chamber
546     nClusCh[ipCh] = nMipClusTot;
547     
548     if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;}                                                                         
549     
550     Int_t index=-1;                                                                             //index of the "best" matching cluster
551     
552     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                               //clusters loop
553       
554       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                    //get the cluster
555 // evaluate qThre
556       if(tsRight){
557        if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {
558          qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts);                                         //
559        } else {                                                                                  // in the past just 1 qthre
560          hvsec = pParam->InHVSector(pClu->Y());                                                  //  per chamber
561          if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts);                                     //
562         } 
563       } else qthre = pParam->QCut();
564             
565 //      
566       if(pClu->Q()<qthre) continue;                                                                      //charge compartible with MIP clusters      
567       isOkQcut = kTRUE;
568 //
569       cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y();                                            //get the LORS coordinates of the cluster
570       Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
571       
572       if(dist<dmin) {
573         dmin = dist;
574         index=iClu;
575         bestHmpCluster=pClu;
576       }
577     } // clusters loop
578
579     if(!isOkQcut) {
580       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
581       delete hmpTrk;hmpTrk=0x0; continue;                                                                     
582     }
583
584     Double_t radius = (pParam->Lors2Mars(ipCh,pParam->SizeAllX()/2,pParam->SizeAllY()/2)).Mag(); 
585     
586     if(!AliTracker::PropagateTrackToBxByBz(hmpTrk,radius,pTrk->GetMass(),1,kFALSE)) {delete hmpTrk;hmpTrk=0x0;continue;}
587               
588     if(!hmpTrk->PropagateTo(bestHmpCluster)) {delete hmpTrk;hmpTrk=0x0;continue;}
589
590     Int_t cluSiz = bestHmpCluster->Size();
591     pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0);  //store mip info in any case 
592     pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                             //set chamber, index of cluster + cluster size
593
594     if(AliHMPIDReconstructor::GetRecoParam())                                                 //retrieve distance cut
595     {
596       if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE)                      //distance cut is fixed number
597       { 
598         distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist();
599       }
600       else 
601       {
602         for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar);      //prevision: distance cut is function of momentum
603         distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise
604       }
605     }
606     else {distCut=pParam->DistCut();}
607       
608     if(dmin < distCut) {
609       isOkDcut = kTRUE;
610     }   
611     
612     if(!isOkDcut) {
613       pTrk->SetHMPIDsignal(pParam->kMipDistCut);                                                //closest cluster with enough charge is still too far from intersection
614     }
615     
616     if(isOkQcut*isOkDcut) isMatched = kTRUE;                                                    // MIP-Track matched !!    
617     
618     if(!isMatched) {delete hmpTrk;hmpTrk=0x0;continue;}                                           // If matched continue...
619
620     Bool_t isOk = hmpTrk->Update(bestHmpCluster,0.1,0);
621     if(!isOk) {delete hmpTrk;hmpTrk=0x0;continue;}
622     pTrk->SetOuterHmpParam(hmpTrk,AliESDtrack::kHMPIDout);                 
623
624     FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
625  
626     Int_t iRad     = pParam->Radiator(yRa);                                                      //evaluate the radiator involved
627     
628     //evaluate nMean
629     if(tsRight){
630      if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
631        nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts);                                               //C6F14 Nmean for this chamber
632      } else {
633        if(iRad < 0) {
634         nmean = -1;
635        } else {
636        Double_t tLow  = ((TF1*)pNmean->At(6*ipCh+2*iRad  ))->Eval(ts);                           //C6F14 low  temp for this chamber
637        Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts);                           //C6F14 high temp for this chamber
638        Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
639        nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
640        }
641        if(nmean < 0){                                                                            //track didn' t pass through the radiator
642          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
643          pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                          //set index of cluster
644          delete hmpTrk;hmpTrk=0x0; 
645          continue;
646         }
647       }
648     } else nmean = pParam->MeanIdxRad();
649     
650     //
651     recon.SetImpPC(xPc,yPc);                                                                     //store track impact to PC
652     recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa);                   //search for Cerenkov angle of this track
653
654     if(pTrk->GetHMPIDsignal()<0) {delete hmpTrk;hmpTrk=0x0;continue;}
655         
656     AliHMPIDPid pID;
657     Double_t prob[5];
658     pID.FindPid(pTrk,nmean,5,prob);
659   //  Printf("Tracker RefIndex = %f, chmaber = %i, prob1 = %f, prob2 = %f, prob3 = %f, prob4 = %f, prob5 = %f", nmean,ipCh,prob[0],prob[1],prob[2],prob[3],prob[4]);
660     pTrk->SetHMPIDpid(prob);
661     delete hmpTrk;hmpTrk=0x0;
662   }//iTrk
663
664   return 0; // error code: 0=no error;
665 }//ReconFromKin()
666 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
667 Int_t AliHMPIDTracker::ReconHiddenTrk(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
668 {
669 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
670 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time), pQthre - pointer to all function Qthre=f(time)
671 //   Returns: error code, 0 if no errors
672   
673   AliHMPIDReconHTA reconHTA;                                                                     //instance of reconstruction class, nothing important in ctor
674   
675   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
676  
677   Bool_t tsRight = kTRUE;
678  
679   UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin();                                        //
680   UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax();                                        //
681   UInt_t ts = pEsd->GetTimeStamp();
682
683   if(ts<tsmin || ts>tsmax) {
684     AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
685     tsRight = kFALSE;
686   }
687    
688   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
689     
690     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //here it is simulated or just empty track
691     Int_t ipCh;
692     ipCh = pTrk->GetHMPIDcluIdx();ipCh/=1000000;
693     if(ipCh<0) continue;
694
695     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
696     Int_t nMipClusTot = pMipCluLst->GetEntries();                                               //total number of clusters in the given chamber
697     
698     Double_t qMip=-1;
699     Int_t chMip=-1;    
700     Double_t xMip = 0;
701     Double_t yMip = 0;
702     Int_t indMip  = 0;
703     Int_t cluMipSiz = 0;
704
705     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                               //clusters loop
706       
707       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                    //get the cluster
708       Double_t qClus = pClu->Q();
709       if(qClus>qMip) {
710         qMip  = qClus;
711         chMip = pClu->Ch();
712         xMip = pClu->X();
713         yMip = pClu->Y();
714         indMip = iClu;
715         cluMipSiz = pClu->Size();
716       }
717     }//clus loop
718
719     if(chMip<0) return 1;    
720     
721     Int_t hvsec;
722     Double_t qthre=0;
723     
724 // evaluate qThre
725     if(tsRight){
726     if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {                                              // just for backward compatibility
727       qthre=((TF1*)pQthre->At(chMip))->Eval(ts);                                                  //
728     } else {                                                                                      // in the past just 1 qthre
729       hvsec = pParam->InHVSector(yMip);                                                           //  per chamber
730       if(hvsec>=0) qthre=((TF1*)pQthre->At(6*chMip+hvsec))->Eval(ts);                             //
731     } 
732    } else qthre = pParam->QCut();
733 //
734     if(qMip<qthre) {
735       pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0);                                                 //store mip info in any case 
736       pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);                                                          
737       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
738       return 1;                                                                                   //charge compatible with MIP clusters      
739     }
740       
741     pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0);                                                   //store mip info in any case 
742     pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);                                                          
743
744     Double_t yRa = yMip;                                                                        //just an approx...
745     Double_t nmean;
746
747     Int_t iRad     = pParam->Radiator(yRa);                                                   //evaluate the radiator involved
748     
749     //evaluate nMean
750     if(tsRight){
751     if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
752       nmean=((TF1*)pNmean->At(3*chMip))->Eval(ts);                                              //C6F14 Nmean for this chamber
753     } else {
754       if(iRad < 0) {
755         nmean = -1;
756       } else {
757       Double_t tLow  = ((TF1*)pNmean->At(6*chMip+2*iRad  ))->Eval(ts);                          //C6F14 low  temp for this chamber
758       Double_t tHigh = ((TF1*)pNmean->At(6*chMip+2*iRad+1))->Eval(ts);                          //C6F14 high temp for this chamber
759       Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
760       nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
761       }
762       if(nmean < 0){                                                                            //track didn' t pass through the radiator
763          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
764          return 1;
765         }
766       } 
767     } else nmean = pParam->MeanIdxRad();
768     //
769     if(!reconHTA.CkovHiddenTrk(pTrk,(TClonesArray *)pClus->At(ipCh),indMip,nmean)) {                 //search for track parameters and Cerenkov angle of this track
770       AliHMPIDPid pID;
771       Double_t prob[5];
772       pID.FindPid(pTrk,nmean,5,prob);
773       pTrk->SetHMPIDpid(prob);
774     }
775 //      Printf(" Prob e- %6.2f mu %6.2f pi %6.2f k %6.2f p %6.2f",prob[0]*100,prob[1]*100,prob[2]*100,prob[3]*100,prob[4]*100);
776   }//iTrk
777
778   return 0; // error code: 0=no error;
779
780 }//Recon()
781 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
782 void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
783   
784  // Publishes all pointers to clusters known to the tracker into the
785   // passed object array.
786   // The ownership is not transfered - the caller is not expected to delete
787   // the clusters
788  
789   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){    
790     TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
791     for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
792       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);    
793       array->AddLast(pClu);
794     }//cluster loop in iCh
795     pCluArr->Delete();
796   }//Ch loop
797     
798   return;
799 }
800 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++