]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDTracker.cxx
Tracking procedure updated, it is based now on AliTrackerBase method. The parameters...
[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   AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
64   Float_t xyz[3];
65   pClu->GetGlobalXYZ(xyz);
66   Float_t cov[6];
67   pClu->GetGlobalCov(cov);
68   point.SetXYZ(xyz,cov);
69   return kTRUE;
70 }//GetTrackPoint()
71 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
72 Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
73 {
74 // Static method to find intersection in between given track and HMPID chambers
75 // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
76 //   Returns: intersected chamber ID or -1
77   AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
78   for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){                              //chambers loop
79     Int_t chInt = IntTrkCha(i,hmpTrk,xPc,yPc,xRa,yRa,theta,phi);
80     if(chInt>=0) {delete hmpTrk;hmpTrk=0x0;return chInt;}
81   }                                                                                             //chambers loop
82   delete hmpTrk; hmpTrk=0x0;
83   return -1;                                                                                    //no intersection with HMPID chambers
84 }//IntTrkCha()
85 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
86 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)
87 {
88 // Static method to find intersection in between given track and HMPID chambers
89 // Arguments: pTrk- HMPID track; xPc,yPc- track intersection with PC in LORS [cm]
90 //   Returns: intersected chamber ID or -1
91     AliHMPIDParam *pParam=AliHMPIDParam::Instance();
92     Double_t p1[3],n1[3]; 
93     pParam->Norm(ch,n1); 
94     pParam->Point(ch,p1,AliHMPIDParam::kRad);                                                    //point & norm  for middle of radiator plane
95     Double_t p2[3],n2[3]; 
96     pParam->Norm(ch,n2);
97     pParam->Point(ch,p2,AliHMPIDParam::kPc);                                                     //point & norm  for entrance to PC plane
98     if(pTrk->Intersect(p1,n1)==kFALSE) return -1;                                                //try to intersect track with the middle of radiator
99     if(pTrk->Intersect(p2,n2)==kFALSE) return -1;   
100     pParam->Mars2LorsVec(ch,n1,theta,phi);                                                       //track angles at RAD
101     pParam->Mars2Lors   (ch,p1,xRa,yRa);                                                         //TRKxRAD position
102     pParam->Mars2Lors   (ch,p2,xPc,yPc);                                                         //TRKxPC position
103     if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return ch;                     //return intersected chamber  
104   return -1;                                                                                     //no intersection with HMPID chambers
105 }//IntTrkCha()
106 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
107 Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
108 {
109 // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
110 // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
111 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
112   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
113   pCluTree->GetEntry(0);
114   return 0;  
115 }//LoadClusters()
116 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
117 Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
118 {
119 // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
120 // Agruments: pEsd - pointer to ESD
121 //   Returns: error code    
122   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
123   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
124   if(!pNmeanEnt) AliError("No Nmean C6F14 ");
125   if(!pQthreEnt) AliError("No Qthre");
126     
127   return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());  
128 }//PropagateBack()
129 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
130 Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
131 {
132 // Static method to reconstruct the Cherenkov angle for all valid tracks of a given event.
133 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
134 // Returns: error code, 0 if no errors   
135 //
136 // Algortihm: Loop over tracks
137 //
138 // 1. Find the closest MIP cluster using fast tracks extrapolation method
139 // 2. Propagate track to the MIP cluster using the STEER method
140 // 3. Update the track information with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
141 // 4. Propagate back the constrained track to the radiator radius ( not exact yet)
142 // 5. Propagation in the last 10 cm with the fast method 
143 // 6. Set ESDtrack information
144 // 7. Calculate the Cherenkov angle
145
146   const Double_t kMaxSnp=0.9;   //maximal snp for prolongation
147   const Double_t kdRadiator=10; // distance between radiator and the plane
148   AliHMPIDRecon recon;                                                                           //instance of reconstruction class, nothing important in ctor
149   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
150   Float_t xPc,yPc,xRa,yRa,theta,phi;
151   
152   Double_t cluLORS[2]={0};
153   Int_t nMipClusTot=0;
154   
155   Double_t qthre = 0;   Double_t nmean=0; Int_t hvsec=0;
156   Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
157
158   Bool_t tsRight = kTRUE;
159   
160   UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin();                                        //
161   UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax();                                        //
162   UInt_t ts = pEsd->GetTimeStamp();
163   
164   if(ts<tsmin || ts>tsmax) {
165     AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
166     tsRight = kFALSE;
167   }
168    
169   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
170            
171     //    Double_t bestChi2=99999;chi2=99999;                                                     //init. track matching params
172     Double_t dmin=999999,bz=0,distCut=1,distParams[5]={1};
173
174     Bool_t isOkDcut=kFALSE;
175     Bool_t isOkQcut=kFALSE;
176     Bool_t isMatched=kFALSE;
177     
178     AliHMPIDCluster *bestHmpCluster=0x0;                                                          //the best matching cluster
179     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //get reconstructed track   
180     
181     if(!pTrk->IsOn(AliESDtrack::kTPCout)) continue;
182  
183     if(pTrk->IsOn(AliESDtrack::kTPCrefit)) continue;
184     AliESDfriendTrack *ftrack= (AliESDfriendTrack *)pTrk->GetFriendTrack();
185     if (!ftrack) continue; 
186     if (!ftrack->GetTPCOut()) continue;
187     AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
188     AliHMPIDtrack *hmpTrkConstrained = 0;                                                         //create a hmpid track to be used for propagation and matching 
189     hmpTrk->Set(ftrack->GetTPCOut()->GetX(), ftrack->GetTPCOut()->GetAlpha(),ftrack->GetTPCOut()->GetParameter(), ftrack->GetTPCOut()->GetCovariance());
190     //
191     bz=AliTracker::GetBz();
192           
193     //initial flags for HMPID ESD infos    
194     pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
195     pTrk->SetHMPIDmip(0,0,0,0);                                                                //store mip info in any case 
196     pTrk->SetHMPIDcluIdx(99,99999);                                                            //chamber not found, mip not yet considered
197     pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
198     
199     Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);                                        //find the intersected chamber for this track 
200     if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;}                                                                         //no intersection at all, go after next track
201
202     pTrk->SetHMPIDtrk(xPc,yPc,theta,phi);                                                        //store initial infos
203     pTrk->SetHMPIDcluIdx(ipCh,9999);                                                             //set chamber, index of cluster + cluster size
204     
205     // track intersects the chamber ipCh: find the MIP          
206     
207     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
208     nMipClusTot = pMipCluLst->GetEntries();                                                     //total number of clusters in the given chamber
209     nClusCh[ipCh] = nMipClusTot;
210     
211     if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;}                                                                         
212     
213     Int_t index=-1;                                                                               //index of the "best" matching cluster
214     //
215     // 1. Find the closest MIP cluster using fast tracks extrapolation method
216     //
217     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                                 //clusters loop
218       
219       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                      //get the cluster
220       // evaluate qThre
221       if(tsRight){
222         if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {
223           qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts);                                          //
224         } else {                                                                                   // in the past just 1 qthre
225           hvsec = pParam->InHVSector(pClu->Y());                                                   //  per chamber
226           if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts);                                     //
227         } 
228       } else qthre = pParam->QCut();
229       
230       //      
231       if(pClu->Q()<qthre) continue;                                                                      //charge compartible with MIP clusters      
232       isOkQcut = kTRUE;
233       //
234       cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y();                                                        //get the LORS coordinates of the cluster
235       Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
236       
237       if(dist<dmin) {
238         dmin = dist;
239         index=iClu;
240         bestHmpCluster=pClu;
241       }
242     } // clusters loop
243
244     // moved down
245     /*if(!isOkQcut) {
246       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
247       delete hmpTrk;hmpTrk=0x0; continue;                                                                     
248     }*/
249     
250     //
251     // 2. Propagate track to the MIP cluster using the STEER method
252     //
253     TVector3 vG = pParam->Lors2Mars(ipCh,bestHmpCluster->X(),bestHmpCluster->Y());    
254     Double_t gx = vG[0];
255     Double_t gy = vG[1];
256     Double_t gz = vG[2];
257     Double_t alpha=TMath::ATan2(gy,gx);
258     Double_t radiusH=TMath::Sqrt(gy*gy+gx*gx);
259     if (!(hmpTrk->Rotate(alpha,kTRUE))) continue;
260     if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrk,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}    
261     //
262     // 3. Update the track with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
263     //
264     AliExternalTrackParam trackC(*hmpTrk);
265     Double_t posC[2]={0,gz};
266     Double_t covC[3]={0.1*0.1, 0, 0.1*0.1};
267     trackC.Update(posC,covC);
268     //
269     // 4. Propagate back the constrained track to the radiator radius ( not exact yet)
270     //
271     hmpTrkConstrained = new AliHMPIDtrack(*pTrk);
272     hmpTrkConstrained->Set(trackC.GetX(), trackC.GetAlpha(),trackC.GetParameter(), trackC.GetCovariance());
273     if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrkConstrained,radiusH-kdRadiator,pTrk->GetMass(),1,kFALSE,kMaxSnp,1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;}
274     //
275     // 5. Propagation in the last 10 cm with the fast method 
276     //
277     Float_t xPc0=0, yPc0=0;
278     
279     IntTrkCha(ipCh, hmpTrk, xPc0,yPc0,xRa,yRa,theta,phi);  
280     IntTrkCha(ipCh, hmpTrkConstrained, xPc,yPc,xRa,yRa,theta,phi);  
281     //
282     // 6. Set ESDtrack information
283     //
284     Int_t cluSiz = bestHmpCluster->Size();
285     pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0);  //store mip info in any case 
286     pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                             //set chamber, index of cluster + cluster size    
287     pTrk->SetHMPIDtrk(xPc0,yPc0,theta,phi);
288     //
289     //
290     // Dump debug info if specified
291     // 
292     if (AliHMPIDReconstructor::StreamLevel()>0) {       
293       AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
294       AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);      
295       trackTPC->Rotate(alpha);
296       trackCurrent->Rotate(alpha);      
297       Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);      
298       Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);    
299       Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
300       Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());    
301       //
302       AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
303       Double_t pos[2]={0,gz};
304       Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
305       Double_t chi2C =  trackTPCConstrained->GetPredictedChi2(pos,cov);
306       trackTPCConstrained->Update(pos,cov);
307       (*fDebugStreamer)<<"track"<<
308         "rH="<<radiusH<<                      // radius of cluster
309         "angle="<<tanAlpha<<                  // tan of the local inlination angle
310         "dC="<<deltaC<<                       // delta of the curvature
311         "trackTPC.="<<trackTPC<<              // TPc outer param extrapolated to the HMPID
312         "chi2C="<<chi2C<<
313         "trackTPCC.="<<trackTPCConstrained<<  // TPc outer param extrapolated to the HMPID
314         "trackCurrent.="<<trackCurrent<<      // current track extrapolated to the HMPID
315         "sTPC="<<statusTPC<<                  // status for the current TPC  track
316         "sCurrent="<<statusCurrent<<          // status for the current global track
317         "cl.="<<bestHmpCluster<<              // HMPID cluster
318         //
319         "t.="<<pTrk<<                        // curent esd track
320         "ft.="<<ftrack<<                     // friend track
321         "hmpTrk.="<<hmpTrk<<                 // hmpid tracks as used in the following code
322         "hmpTrkC.="<<hmpTrkConstrained<<     // hmpid tracks as used in the following code
323         "gx="<<gx<<                          // global cluster position X
324         "gy="<<gy<<                          // Y
325         "gz="<<gz<<                          // Z
326         "\n";
327     }                 
328     //
329     //
330     //    
331     if(!isOkQcut) {
332       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
333       delete hmpTrk;hmpTrk=0x0; 
334       delete hmpTrkConstrained;hmpTrkConstrained=0x0; 
335       continue;                                                                     
336     }        
337     
338     if(AliHMPIDReconstructor::GetRecoParam())                                                 //retrieve distance cut
339     {
340       if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE)                      //distance cut is fixed number
341       { 
342         distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist();
343       }
344       else 
345       {
346         for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar);      //prevision: distance cut is function of momentum
347         distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise
348       }
349     }
350     else {distCut=pParam->DistCut();}
351      
352     //dmin recalculated
353     
354     dmin = TMath::Sqrt((xPc0-bestHmpCluster->X())*(xPc0-bestHmpCluster->X())+(yPc0-bestHmpCluster->Y())*(yPc0-bestHmpCluster->Y()));
355              
356     if(dmin < distCut) {
357       isOkDcut = kTRUE;
358     }   
359     //isOkDcut = kTRUE; // switch OFF cut
360     
361     if(!isOkDcut) {
362       pTrk->SetHMPIDsignal(pParam->kMipDistCut);                                                //closest cluster with enough charge is still too far from intersection
363     }
364     
365     if(isOkQcut*isOkDcut) isMatched = kTRUE;                                                    // MIP-Track matched !!    
366     
367     if(!isMatched) {delete hmpTrk;hmpTrk=0x0;delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;}                                           // If matched continue...
368
369     Bool_t isOk = kTRUE; 
370     if(!isOk) {delete hmpTrk;hmpTrk=0x0;continue; delete hmpTrkConstrained;hmpTrkConstrained=0x0;}
371     pTrk->SetOuterHmpParam(hmpTrkConstrained,AliESDtrack::kHMPIDout);                 
372
373     FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
374  
375     Int_t iRad = pParam->Radiator(yRa);                                                          //evaluate the radiator involved
376     
377     //evaluate nMean
378     if(tsRight){
379      if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
380        nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts);                                               //C6F14 Nmean for this chamber
381      } else {
382        if(iRad < 0) {
383         nmean = -1;
384        } else {
385        Double_t tLow  = ((TF1*)pNmean->At(6*ipCh+2*iRad  ))->Eval(ts);                           //C6F14 low  temp for this chamber
386        Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts);                           //C6F14 high temp for this chamber
387        Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
388        nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
389        }
390        if(nmean < 0){                                                                            //track didn' t pass through the radiator
391          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
392          pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                          //set index of cluster
393          delete hmpTrk;hmpTrk=0x0; 
394          delete hmpTrkConstrained;hmpTrkConstrained=0x0; 
395          continue;
396         }
397      }
398     } else nmean = pParam->MeanIdxRad();    
399     //
400     // 7. Calculate the Cherenkov angle
401     //
402     recon.SetImpPC(xPc0,yPc0);                                                                     //store track impact to PC           
403     recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa);                   //search for Cerenkov angle of this track
404     
405     Double_t thetaCkov = pTrk->GetHMPIDsignal();
406     
407     if (AliHMPIDReconstructor::StreamLevel()>0) {       
408       AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
409       AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);      
410       trackTPC->Rotate(alpha);
411       trackCurrent->Rotate(alpha);      
412       Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);      
413       Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);    
414       Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
415       Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());    
416       //
417       AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
418       Double_t pos[2]={0,gz};
419       Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
420       Double_t chi2C =  trackTPCConstrained->GetPredictedChi2(pos,cov);
421       trackTPCConstrained->Update(pos,cov);
422       (*fDebugStreamer)<<"track2"<<
423         "rH="<<radiusH<<                      // radius of cluster
424         "angle="<<tanAlpha<<                  // tan of the local inlination angle
425         "dC="<<deltaC<<                       // delta of the curvature
426         "trackTPC.="<<trackTPC<<              // TPc outer param extrapolated to the HMPID
427         "chi2C="<<chi2C<<
428         "trackTPCC.="<<trackTPCConstrained<<  // TPc outer param extrapolated to the HMPID
429         "trackCurrent.="<<trackCurrent<<      // current track extrapolated to the HMPID
430         "sTPC="<<statusTPC<<                  // status for the current TPC  track
431         "sCurrent="<<statusCurrent<<          // status for the current global track
432         "cl.="<<bestHmpCluster<<              // HMPID cluster
433         //
434         "t.="<<pTrk<<                         // curent esd track
435         "ft.="<<ftrack<<                      // friend track
436         "hmpTrk.="<<hmpTrk<<                  // hmpid tracks as used in the following code
437         "hmpTrkC.="<<hmpTrkConstrained<<      // hmpid tracks as used in the following code
438         "gx="<<gx<<                           // global cluster position X
439         "gy="<<gy<<                           // Y
440         "gz="<<gz<<                           // Z
441         "thetaCkov="<<thetaCkov<<             // Cherenkov angle
442         "\n";
443     } 
444                 
445     if(pTrk->GetHMPIDsignal()<0) {
446       delete hmpTrk;hmpTrk=0x0;
447       delete hmpTrkConstrained; hmpTrkConstrained=0x0;
448       continue;}
449         
450     AliHMPIDPid pID;
451     Double_t prob[5];
452     pID.FindPid(pTrk,5,prob);
453     pTrk->SetHMPIDpid(prob);
454 //      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);
455     delete hmpTrk; hmpTrk=0x0;
456     delete hmpTrkConstrained; hmpTrkConstrained=0x0;
457   }//iTrk
458
459   return 0; // error code: 0=no error;
460 }//Recon()
461 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
462 Int_t AliHMPIDTracker::ReconHiddenTrk(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
463 {
464 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
465 // 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)
466 //   Returns: error code, 0 if no errors
467   
468   AliHMPIDReconHTA reconHTA;                                                                     //instance of reconstruction class, nothing important in ctor
469   
470   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
471  
472   Bool_t tsRight = kTRUE;
473  
474   UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin();                                        //
475   UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax();                                        //
476   UInt_t ts = pEsd->GetTimeStamp();
477
478   if(ts<tsmin || ts>tsmax) {
479     AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
480     tsRight = kFALSE;
481   }
482    
483   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
484     
485     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //here it is simulated or just empty track
486     Int_t ipCh;
487     ipCh = pTrk->GetHMPIDcluIdx();ipCh/=1000000;
488     if(ipCh<0) continue;
489
490     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
491     Int_t nMipClusTot = pMipCluLst->GetEntries();                                               //total number of clusters in the given chamber
492     
493     Double_t qMip=-1;
494     Int_t chMip=-1;    
495     Double_t xMip = 0;
496     Double_t yMip = 0;
497     Int_t indMip  = 0;
498     Int_t cluMipSiz = 0;
499
500     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                               //clusters loop
501       
502       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                    //get the cluster
503       Double_t qClus = pClu->Q();
504       if(qClus>qMip) {
505         qMip  = qClus;
506         chMip = pClu->Ch();
507         xMip = pClu->X();
508         yMip = pClu->Y();
509         indMip = iClu;
510         cluMipSiz = pClu->Size();
511       }
512     }//clus loop
513
514     if(chMip<0) return 1;    
515     
516     Int_t hvsec;
517     Double_t qthre=0;
518     
519 // evaluate qThre
520     if(tsRight){
521     if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {                                              // just for backward compatibility
522       qthre=((TF1*)pQthre->At(chMip))->Eval(ts);                                                  //
523     } else {                                                                                      // in the past just 1 qthre
524       hvsec = pParam->InHVSector(yMip);                                                           //  per chamber
525       if(hvsec>=0) qthre=((TF1*)pQthre->At(6*chMip+hvsec))->Eval(ts);                             //
526     } 
527    } else qthre = pParam->QCut();
528 //
529     if(qMip<qthre) {
530       pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0);                                                 //store mip info in any case 
531       pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);                                                          
532       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
533       return 1;                                                                                   //charge compatible with MIP clusters      
534     }
535       
536     pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0);                                                   //store mip info in any case 
537     pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);                                                          
538
539     Double_t yRa = yMip;                                                                        //just an approx...
540     Double_t nmean;
541
542     Int_t iRad     = pParam->Radiator(yRa);                                                   //evaluate the radiator involved
543     
544     //evaluate nMean
545     if(tsRight){
546     if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
547       nmean=((TF1*)pNmean->At(3*chMip))->Eval(ts);                                              //C6F14 Nmean for this chamber
548     } else {
549       if(iRad < 0) {
550         nmean = -1;
551       } else {
552       Double_t tLow  = ((TF1*)pNmean->At(6*chMip+2*iRad  ))->Eval(ts);                          //C6F14 low  temp for this chamber
553       Double_t tHigh = ((TF1*)pNmean->At(6*chMip+2*iRad+1))->Eval(ts);                          //C6F14 high temp for this chamber
554       Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
555       nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
556       }
557       if(nmean < 0){                                                                            //track didn' t pass through the radiator
558          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
559          return 1;
560         }
561       } 
562     } else nmean = pParam->MeanIdxRad();
563     //
564     if(!reconHTA.CkovHiddenTrk(pTrk,(TClonesArray *)pClus->At(ipCh),indMip,nmean)) {                 //search for track parameters and Cerenkov angle of this track
565       AliHMPIDPid pID;
566       Double_t prob[5];
567       pID.FindPid(pTrk,5,prob);
568       pTrk->SetHMPIDpid(prob);
569     }
570 //      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);
571   }//iTrk
572
573   return 0; // error code: 0=no error;
574
575 }//Recon()
576 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
577 void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
578   
579  // Publishes all pointers to clusters known to the tracker into the
580   // passed object array.
581   // The ownership is not transfered - the caller is not expected to delete
582   // the clusters
583  
584   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){    
585     TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
586     for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
587       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);    
588       array->AddLast(pClu);
589     }//cluster loop in iCh
590     pCluArr->Delete();
591   }//Ch loop
592     
593   return;
594 }
595 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++