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