]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDTracker.cxx
e7dcbb02ed184e711b54c0073598164315e18a7f
[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 "AliHMPIDReconHTA.h"    //ReconHTA()
8 #include <AliESDEvent.h>         //PropagateBack(),Recon()  
9 #include <AliESDtrack.h>         //Intersect() 
10 #include <AliTracker.h> 
11 #include <AliRun.h>              //GetTrackPoint(),PropagateBack()  
12 #include <AliTrackPointArray.h>  //GetTrackPoint()
13 #include <AliAlignObj.h>         //GetTrackPoint()
14 #include <AliCDBManager.h>       //PropageteBack()
15 #include <AliCDBEntry.h>         //PropageteBack()
16 //.
17 // HMPID base class fo tracking
18 //.
19 //.
20 //.
21 ClassImp(AliHMPIDTracker)
22 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 AliHMPIDTracker::AliHMPIDTracker():
24   AliTracker(),
25   fClu(new TObjArray(AliHMPIDParam::kMaxCh+1))  
26 {
27 // ctor. Create TObjArray of TClonesArray of AliHMPIDCluster  
28 // 
29 //  
30   fClu->SetOwner(kTRUE);
31   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
32 }//ctor
33 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++    
34 Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
35 {
36 // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
37 // MIP cluster is reffered by index which is stored in AliESDtrack  ???????
38 // Arguments: idx- cluster index which is stored by HMPID in AliESDtrack
39 //            point- reference to the object where to store the point     
40 //   Returns: status of operation  if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. 
41   if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
42   Int_t iCham=idx/1000000; Int_t iClu=idx%1000000;
43   iClu = iClu%1000; //GetHMPIDcluIdx -> 1e+6*ch + 1e+3*clusize + cluIdx;
44   point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham));//layer and chamber number
45   TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham];
46   AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
47   Float_t xyz[3];
48   pClu->GetGlobalXYZ(xyz);
49   Float_t cov[6];
50   pClu->GetGlobalCov(cov);
51   point.SetXYZ(xyz,cov);
52   return kTRUE;
53 }//GetTrackPoint()
54 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
55 Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
56 {
57 // Static method to find intersection in between given track and HMPID chambers
58 // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
59 //   Returns: intersected chamber ID or -1
60   AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
61   for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){                              //chambers loop
62     Int_t chInt = IntTrkCha(i,hmpTrk,xPc,yPc,xRa,yRa,theta,phi);
63     if(chInt>=0) {delete hmpTrk;hmpTrk=0x0;return chInt;}
64   }                                                                                             //chambers loop
65   delete hmpTrk; hmpTrk=0x0;
66   return -1;                                                                                    //no intersection with HMPID chambers
67 }//IntTrkCha()
68 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
69 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)
70 {
71 // Static method to find intersection in between given track and HMPID chambers
72 // Arguments: pTrk- HMPID track; xPc,yPc- track intersection with PC in LORS [cm]
73 //   Returns: intersected chamber ID or -1
74     AliHMPIDParam *pParam=AliHMPIDParam::Instance();
75     Double_t p1[3],n1[3]; 
76     pParam->Norm(ch,n1); 
77     pParam->Point(ch,p1,AliHMPIDParam::kRad);                                                    //point & norm  for middle of radiator plane
78     Double_t p2[3],n2[3]; 
79     pParam->Norm(ch,n2);
80     pParam->Point(ch,p2,AliHMPIDParam::kPc);                                                     //point & norm  for entrance to PC plane
81     if(pTrk->Intersect(p1,n1)==kFALSE) return -1;                                                //try to intersect track with the middle of radiator
82     if(pTrk->Intersect(p2,n2)==kFALSE) return -1;   
83     pParam->Mars2LorsVec(ch,n1,theta,phi);                                                       //track angles at RAD
84     pParam->Mars2Lors   (ch,p1,xRa,yRa);                                                         //TRKxRAD position
85     pParam->Mars2Lors   (ch,p2,xPc,yPc);                                                         //TRKxPC position
86     if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return ch;                     //return intersected chamber  
87   return -1;                                                                                     //no intersection with HMPID chambers
88 }//IntTrkCha()
89 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
90 Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
91 {
92 // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
93 // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
94 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
95   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
96   pCluTree->GetEntry(0);
97   return 0;  
98 }//LoadClusters()
99 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
100 Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
101 {
102 // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
103 // Agruments: pEsd - pointer to ESD
104 //   Returns: error code    
105   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
106   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
107   if(!pNmeanEnt) AliError("No Nmean C6F14 ");
108   if(!pQthreEnt) AliError("No Qthre");
109     
110   return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());  
111 }//PropagateBack()
112 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113 Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
114 {
115 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
116 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
117 //   Returns: error code, 0 if no errors   
118   
119   AliHMPIDRecon recon;                                                                           //instance of reconstruction class, nothing important in ctor
120   Float_t xPc,yPc,xRa,yRa,theta,phi;
121   Double_t cluLORS[2]={0};
122 //  Double_t cluMARS[3]={0},trkMARS[3]={0};
123 //  Double_t bestcluMARS[3]={0,0,0};
124 //  Double_t radClu,radInitTrk;   
125   Int_t nMipClusTot=0;
126 //  Double_t d3d=0;
127   Double_t qthre = 0;   Double_t nmean=0; Int_t hvsec=0;
128   Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
129   
130   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
131   
132   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
133            
134 //    Double_t bestChi2=99999;chi2=99999;                                                          //init. track matching params
135     Double_t dmin=999999,bz=0;
136
137     Bool_t isOkDcut=kFALSE;
138     Bool_t isOkQcut=kFALSE;
139     Bool_t isMatched=kFALSE;
140     
141     AliHMPIDCluster *bestHmpCluster=0x0;                                                          //the best matching cluster
142     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //get reconstructed track   
143     
144     if(!pTrk->IsOn(AliESDtrack::kTPCout)) continue;
145  
146     if(pTrk->IsOn(AliESDtrack::kTPCrefit)) continue;
147  
148     AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
149     bz=AliTracker::GetBz();  
150 //initial flags for HMPID ESD infos    
151     pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
152     pTrk->SetHMPIDmip(0,0,0,0);                                                                //store mip info in any case 
153     pTrk->SetHMPIDcluIdx(99,99999);                                                            //chamber not found, mip not yet considered
154     pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
155     
156     Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);                                        //find the intersected chamber for this track 
157
158     if(ipCh<0) continue;                                                                         //no intersection at all, go after next track
159
160     pTrk->SetHMPIDtrk(xPc,yPc,theta,phi);                                                        //store initial infos
161     pTrk->SetHMPIDcluIdx(ipCh,9999);                                                             //set chamber, index of cluster + cluster size
162     
163 // track intersects the chamber ipCh: find the MIP          
164     
165     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
166     nMipClusTot = pMipCluLst->GetEntries();                                                     //total number of clusters in the given chamber
167     nClusCh[ipCh] = nMipClusTot;
168     
169     if(nMipClusTot==0) continue;                                                                         
170     
171     Int_t index=-1;                                                                             //index of the "best" matching cluster
172     
173     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                               //clusters loop
174       
175       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                    //get the cluster
176 // evaluate qThre
177       if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {                                             // just for backward compatibility
178         qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(pEsd->GetTimeStamp());                          //
179       } else {                                                                                     // in the past just 1 qthre
180         hvsec = pParam->InHVSector(pClu->Y());                                              //  per chamber
181         if(hvsec>=0)
182           qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(pEsd->GetTimeStamp());                      //
183       }                                                                                            //
184 //
185       if(pClu->Q()<qthre) continue;                                                                      //charge compartible with MIP clusters      
186       isOkQcut = kTRUE;
187 //
188       cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y();                                            //get the LORS coordinates of the cluster
189       Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
190       
191       if(dist<dmin) {
192         dmin = dist;
193         index=iClu;
194         bestHmpCluster=pClu;
195       }
196     }
197 /*            
198       pParam->Lors2Mars(ipCh,cluLORS[0],cluLORS[1],cluMARS);              //convert cluster coors. from LORS to MARS
199       radClu=TMath::Sqrt(cluMARS[0]*cluMARS[0]+cluMARS[1]*cluMARS[1]);                       //radial distance of candidate cluster in MARS                                          
200       Double_t trkx0[3]; 
201       hmpTrk->GetXYZ(trkx0);                                                                 //get track position in MARS
202       radInitTrk=TMath::Sqrt(trkx0[0]*trkx0[0]+trkx0[1]*trkx0[1]);
203       hmpTrk->PropagateToR(radClu,10);
204       hmpTrk->GetXYZ(trkx0);                                                                   //get track position in MARS
205       hmpTrk->GetXYZAt(radClu,bz,trkMARS);                                                     //get the track coordinates at the rad distance after prop. 
206       d3d=TMath::Sqrt((cluMARS[0]-trkMARS[0])*(cluMARS[0]-trkMARS[0])+(cluMARS[1]-trkMARS[1])*(cluMARS[1]-trkMARS[1])+(cluMARS[2]-trkMARS[2])*(cluMARS[2]-trkMARS[2]));
207       chi2=hmpTrk->GetPredictedChi2(pClu);
208       if(dmin > d3d ) {   +                                                                      //to be saved for the moment...
209         cluSiz = pClu->Size();
210         dmin=d3d;
211         bestHmpCluster=pClu;
212         index=iClu;
213         bestChi2=chi2;
214         cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y();
215 //        pParam->Lors2Mars(ipCh,cluLORS[0],cluLORS[1],bestcluMARS); 
216       }//global dmin cut 
217     }//clus loop
218 */
219     if(!isOkQcut) {
220       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
221       continue;                                                                     
222     }
223
224     Double_t radius = (pParam->Lors2Mars(ipCh,pParam->SizeAllX()/2,pParam->SizeAllY()/2)).Mag(); 
225     
226     if(!AliTracker::PropagateTrackTo(hmpTrk,radius,pTrk->GetMass(),1,kFALSE)) continue;
227               
228     if(!hmpTrk->PropagateTo(bestHmpCluster)) continue;
229
230     Int_t cluSiz = bestHmpCluster->Size();
231     pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0);  //store mip info in any case 
232     pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                             //set chamber, index of cluster + cluster size
233     
234     if(dmin < pParam->DistCut()) {
235       isOkDcut = kTRUE;
236     }
237
238     if(!isOkDcut) {
239       pTrk->SetHMPIDsignal(pParam->kMipDistCut);                                                //closest cluster with enough charge is still too far from intersection
240     }
241     
242     if(isOkQcut*isOkDcut) isMatched = kTRUE;                                                    // MIP-Track matched !!    
243     
244     if(!isMatched) continue;                                                                    // If matched continue...
245
246     Bool_t isOk = hmpTrk->Update(bestHmpCluster,0.1,0);
247     if(!isOk) continue;
248     pTrk->SetOuterHmpParam(hmpTrk,AliESDtrack::kHMPIDout);                 
249
250     FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
251     
252     /*    
253     Int_t indexAll = 0;
254     for(Int_t iC=0;iC<ipCh;iC++) indexAll+=nClusCh[iC]; indexAll+=index;                        //to be verified...
255
256     Bool_t isOk = hmpTrk->Update(bestHmpCluster,bestChi2,indexAll);
257     if(!isOk) continue;
258     pTrk->SetOuterParam(hmpTrk,AliESDtrack::kHMPIDout);                 
259
260     cham=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);
261     
262     if(cham<0) {                                                                                  //no intersection at all, go after next track
263       pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
264       pTrk->SetHMPIDcluIdx   (99,99999);                                                         //chamber not found, mip not yet considered
265       pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
266       continue;                                                                         
267     }
268 */
269     //evaluate nMean
270     if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
271       nmean=((TF1*)pNmean->At(3*ipCh))->Eval(pEsd->GetTimeStamp());                             //C6F14 Nmean for this chamber
272     } else {
273       Int_t iRad     = pParam->Radiator(yRa);                                                   //evaluate the radiator involved
274       if(iRad < 0) {
275         nmean = -1;
276       } else {
277       Double_t tLow  = ((TF1*)pNmean->At(6*ipCh+2*iRad  ))->Eval(pEsd->GetTimeStamp());         //C6F14 low  temp for this chamber
278       Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(pEsd->GetTimeStamp());         //C6F14 high temp for this chamber
279       Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
280       nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
281       }
282       if(nmean < 0){                                                                            //track didn' t pass through the radiator
283          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
284          pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                          //set index of cluster
285          continue;
286       }
287     }
288     //
289     recon.SetImpPC(xPc,yPc);                                                                     //store track impact to PC
290     recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa);                   //search for Cerenkov angle of this track
291     
292     AliHMPIDPid pID;
293     Double_t prob[5];
294     pID.FindPid(pTrk,5,prob);
295     pTrk->SetHMPIDpid(prob);
296 //      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);
297
298     delete hmpTrk;hmpTrk=0x0;
299   }//iTrk
300
301   return 0; // error code: 0=no error;
302 }//Recon()
303 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
304 Int_t AliHMPIDTracker::ReconHiddenTrk(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
305 {
306 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
307 // 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)
308 //   Returns: error code, 0 if no errors
309   
310   AliHMPIDReconHTA reconHTA;                                                                     //instance of reconstruction class, nothing important in ctor
311   
312   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
313   
314   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
315     
316     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //here it is simulated or just empty track
317     Int_t ipCh;
318     ipCh = pTrk->GetHMPIDcluIdx();ipCh/=1000000;
319     if(ipCh<0) continue;
320
321     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
322     Int_t nMipClusTot = pMipCluLst->GetEntries();                                               //total number of clusters in the given chamber
323     
324     Double_t qMip=-1;
325     Int_t chMip=-1;    
326     Double_t xMip = 0;
327     Double_t yMip = 0;
328     Int_t indMip  = 0;
329     Int_t cluMipSiz = 0;
330
331     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                               //clusters loop
332       
333       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                    //get the cluster
334       Double_t qClus = pClu->Q();
335       if(qClus>qMip) {
336         qMip  = qClus;
337         chMip = pClu->Ch();
338         xMip = pClu->X();
339         yMip = pClu->Y();
340         indMip = iClu;
341         cluMipSiz = pClu->Size();
342       }
343     }//clus loop
344
345     if(chMip<0) return 1;    
346     
347     Int_t hvsec;
348     Double_t qthre=0;
349 // evaluate qThre
350     if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {                                              // just for backward compatibility
351       qthre=((TF1*)pQthre->At(chMip))->Eval(pEsd->GetTimeStamp());                                //
352     } else {                                                                                      // in the past just 1 qthre
353       hvsec = pParam->InHVSector(yMip);                                                           //  per chamber
354       if(hvsec>=0) qthre=((TF1*)pQthre->At(6*chMip+hvsec))->Eval(pEsd->GetTimeStamp());           //
355     }
356 //
357     if(qMip<qthre) {
358       pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0);                                                        //store mip info in any case 
359       pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);                                                          
360       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
361       return 1;                                                                                   //charge compatible with MIP clusters      
362     }
363       
364     pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0);                                                   //store mip info in any case 
365     pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);                                                          
366
367     Double_t yRa = yMip;                                                                        //just an approx...
368     Double_t nmean;
369     //evaluate nMean
370     if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
371       nmean=((TF1*)pNmean->At(3*chMip))->Eval(pEsd->GetTimeStamp());                            //C6F14 Nmean for this chamber
372     } else {
373       Int_t iRad     = pParam->Radiator(yRa);                                                   //evaluate the radiator involved
374       if(iRad < 0) {
375         nmean = -1;
376       } else {
377       Double_t tLow  = ((TF1*)pNmean->At(6*chMip+2*iRad  ))->Eval(pEsd->GetTimeStamp());        //C6F14 low  temp for this chamber
378       Double_t tHigh = ((TF1*)pNmean->At(6*chMip+2*iRad+1))->Eval(pEsd->GetTimeStamp());        //C6F14 high temp for this chamber
379       Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
380       nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
381       }
382       if(nmean < 0){                                                                            //track didn' t pass through the radiator
383          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
384          return 1;
385       }
386     }
387     //
388     if(!reconHTA.CkovHiddenTrk(pTrk,(TClonesArray *)pClus->At(ipCh),indMip,nmean)) {                 //search for track parameters and Cerenkov angle of this track
389       AliHMPIDPid pID;
390       Double_t prob[5];
391       pID.FindPid(pTrk,5,prob);
392       pTrk->SetHMPIDpid(prob);
393     }
394 //      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);
395   }//iTrk
396
397   return 0; // error code: 0=no error;
398
399 }//Recon()
400 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
401 void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
402   
403  // Publishes all pointers to clusters known to the tracker into the
404   // passed object array.
405   // The ownership is not transfered - the caller is not expected to delete
406   // the clusters
407  
408   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){    
409     TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
410     for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
411       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);    
412       array->AddLast(pClu);
413     }//cluster loop in iCh
414     pCluArr->Delete();
415   }//Ch loop
416     
417   return;
418 }
419 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++