]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDReconHTA.cxx
Merging THbtp and HBTP in one library. Comiplation on Windows/Cygwin
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDReconHTA.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 //////////////////////////////////////////////////////////////////////////
17 //                                                                      //
18 // AliHMPIDReconHTA                                                         //
19 //                                                                      //
20 // HMPID class to perfom pattern recognition based on Hough transfrom    //
21 // for single chamber                                                   //
22 //////////////////////////////////////////////////////////////////////////
23
24 #include "AliHMPIDReconHTA.h"//class header
25 #include "AliHMPIDCluster.h" //CkovHiddenTrk()
26 #include "AliHMPIDRecon.h"   //FunMinPhot()
27 #include <TMinuit.h>         //FitFree()
28 #include <TClonesArray.h>    //CkovHiddenTrk()
29 #include <AliESDtrack.h>     //CkovHiddenTrk()
30 #include <TH2I.h>            //InitDatabase()
31 #include <TGraph.h>          //ShapeModel()
32 #include <TSpline.h>         //ShapeModel()
33 #include "TStopwatch.h"      //
34
35 TH2I* AliHMPIDReconHTA::fgDatabase = 0x0;
36
37 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38 AliHMPIDReconHTA::AliHMPIDReconHTA():
39   TTask("RichRec","RichPat"),
40   fMipX(-999),
41   fMipY(-999),
42   fMipQ(-999),
43   fRadX(-999),
44   fRadY(-999),
45   fIdxMip(0),
46   fNClu(0),
47   fXClu(0),
48   fYClu(0),
49   fClCk(0),
50   fThTrkFit(-999),
51   fPhTrkFit(-999),
52   fCkovFit(-999),
53   fCkovSig2(0),
54   fParam(AliHMPIDParam::Instance())
55 {
56 //..
57 //hidden algorithm
58 //..
59   fParam->SetRefIdx(fParam->MeanIdxRad()); // initialization of ref index to a default one
60   if(!fgDatabase) InitDatabase();
61 }
62 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
63 AliHMPIDReconHTA::~AliHMPIDReconHTA()
64 {
65   DeleteVars();
66 }
67 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
68 void AliHMPIDReconHTA::InitVars(Int_t n)
69 {
70 //..
71 //Init some variables
72 //..
73   fXClu = new Double_t[n];
74   fYClu = new Double_t[n];
75   fClCk = new Bool_t[n];
76   for(Int_t i=0;i<n;i++) fClCk[i] = kTRUE;
77 //
78 }
79 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
80 void AliHMPIDReconHTA::DeleteVars()const
81 {
82 //..
83 //Delete variables
84 //..
85   if(fXClu) delete fXClu;
86   if(fYClu) delete fYClu;
87   if(fClCk) delete fClCk;
88 }
89 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
90 Bool_t AliHMPIDReconHTA::CkovHiddenTrk(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t nmean, Double_t qthre)
91 {
92 // Pattern recognition method without any infos from tracking:HTA (Hidden Track Algorithm)...
93 // The method finds in the chmuber the cluster with the highest charge
94 // compatibile with a MIP, then the strategy is applied
95 // Arguments:  pTrk     - pointer to ESD track
96 //             pCluLs   - list of clusters for a given chamber 
97 //             nmean    - mean freon ref. index
98 //   Returns:           - 0=ok,1=not fitted 
99   
100   AliHMPIDParam *pParam = AliHMPIDParam::Instance(); 
101   pParam->SetRefIdx(nmean);
102
103   if(!CluPreFilter(pCluLst)) {return kFALSE;}
104   
105   Float_t mipX=-1,mipY=-1;Int_t mipId=-1,mipQ=-1;                                                                           
106   Double_t qRef = 0;
107   Int_t nCh=0;
108   for (Int_t iClu=0;iClu<pCluLst->GetEntriesFast();iClu++){                                   //clusters loop
109     AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluLst->UncheckedAt(iClu);                       //get pointer to current cluster    
110     nCh = pClu->Ch();
111     fXClu[iClu] = pClu->X();fYClu[iClu] = pClu->Y();                                          //store x,y for fitting procedure
112     fClCk[iClu] = kTRUE;                                                                      //all cluster are accepted at this stage to be reconstructed
113     if(pClu->Q()>qthre) fClCk[iClu] = kFALSE;                                                 // not a good photon in any case (multiple MIPs)
114     if(pClu->Q()>qRef){                                                                       //searching the highest charge to select a MIP      
115       qRef = pClu->Q();
116       mipId=iClu; mipX=pClu->X();mipY=pClu->Y();mipQ=(Int_t)pClu->Q();
117     }                                                                                    
118 //    Printf(" n. %d x %f y %f Q %f",iClu,pClu->X(),pClu->Y(),pClu->Q());
119   }//clusters loop
120
121   fNClu = pCluLst->GetEntriesFast();
122   if(qRef>qthre){                                                                     //charge compartible with MIP clusters
123     fIdxMip = mipId;
124     fClCk[mipId] = kFALSE;
125     fMipX = mipX; fMipY=mipY; fMipQ = qRef;
126 //    Printf(" mipId %d x %f y %f Q %f",fIdxMip,fMipX,fMipY,fMipQ);
127     if(!DoRecHiddenTrk()) {
128       pTrk->SetHMPIDsignal(pParam->kNoPhotAccept);
129       return kFALSE;
130     }                                                                           //Do track and ring reconstruction,if problems returns 1
131     pTrk->SetHMPIDtrk(fRadX,fRadY,fThTrkFit,fPhTrkFit);                                        //store track intersection info
132     pTrk->SetHMPIDmip(fMipX,fMipY,(Int_t)fMipQ,fNClu);                                         //store mip info 
133     pTrk->SetHMPIDcluIdx(nCh,fIdxMip);                                                         //set cham number and index of cluster
134     pTrk->SetHMPIDsignal(fCkovFit);                                                            //find best Theta ckov for ring i.e. track
135     pTrk->SetHMPIDchi2(fCkovSig2);                                                             //errors squared
136 //    Printf(" n clusters tot %i accepted %i",pCluLst->GetEntriesFast(),fNClu);
137 //    Printf("CkovHiddenTrk: thetaC %f th %f ph %f",fCkovFit,fThTrkFit,fPhTrkFit);
138     return kTRUE;
139   }
140   
141   return kFALSE;
142 }//CkovHiddenTrk()
143 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
144 Bool_t AliHMPIDReconHTA::DoRecHiddenTrk()
145 {
146 // Pattern recognition method without any infos from tracking...
147 // First a preclustering filter to avoid part of the noise
148 // Then only ellipsed-rings are fitted (no possibility, 
149 // for the moment, to reconstruct very inclined tracks)
150 // Finally a fitting with (th,ph) free, starting by very close values
151 // previously evaluated.
152 // Arguments:   none
153 //   Returns:   none
154   Double_t thTrkRec,phiTrkRec,thetaCRec;
155   
156   if(!FindShape(thTrkRec,phiTrkRec,thetaCRec)) return kFALSE;
157
158   if(!FitFree(thTrkRec,phiTrkRec)) {return kFALSE;}
159   
160   return kTRUE;
161 }//DoRecHiddenTrk()
162 /*
163 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164 Bool_t AliHMPIDReconHTA::DoRecHiddenTrk(TClonesArray *pCluLst)
165 {
166 // Pattern recognition method without any infos from tracking...
167 // First a preclustering filter to avoid part of the noise
168 // Then only ellipsed-rings are fitted (no possibility, 
169 // for the moment, to reconstruct very inclined tracks)
170 // Finally a fitting with (th,ph) free, starting by very close values
171 // previously evaluated.
172 // Arguments:   none
173 //   Returns:   none
174   Double_t thTrkRec,phiTrkRec,thetaCRec;
175   
176   if(!FindShape(thTrkRec,phiTrkRec,thetaCRec)) return kFALSE;
177
178 //  Printf("thTrkRec %f phiTrkRec %f ThetaCRec %f",thTrkRec*TMath::RadToDeg(),phiTrkRec*TMath::RadToDeg(),thetaCRec*TMath::RadToDeg());  
179   Int_t nClTmp1 = pCluLst->GetEntriesFast()-1;  //minus MIP...
180   Int_t nClTmp2 = 0;
181   
182   while(nClTmp1 != nClTmp2){
183     SetNClu(pCluLst->GetEntriesFast());
184     if(!FitFree(thTrkRec,phiTrkRec)) {return kFALSE;}
185     nClTmp2 = NClu();
186     if(nClTmp2!=nClTmp1) {nClTmp1=nClTmp2;nClTmp2=0;}
187   }
188   
189   fNClu = nClTmp2;
190   return kTRUE;
191 }//DoRecHiddenTrk()
192 */
193 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
194 Bool_t AliHMPIDReconHTA::CluPreFilter(TClonesArray *pCluLst)
195 {
196 // Pre-filter of bkg clusters
197 // Arguments:    pSluLst  -  List of the clusters for a given chamber
198 //   Returns:    status   -  TRUE if filtering leaves enough photons, FALSE if not
199 //
200   Int_t nClusTot = pCluLst->GetEntriesFast();
201   if(nClusTot<4||nClusTot>100) {
202     return kFALSE; 
203   } else { 
204     InitVars(nClusTot);
205     return kTRUE;
206   }
207 }
208 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
209 Bool_t AliHMPIDReconHTA::FindShape(Double_t &thTrkRec,Double_t &phiTrkRec,Double_t &thetaCRec)
210 {
211 // Finds the estimates for phi and theta of the track and the ThetaCerenkov
212 // by using a database of the shapes of the rings
213 // Arguments:   none  
214 //   Returns:   thTrkRec  - estimate of theta track
215 //              phiTrkRec - estimate of phi   track
216 //              thetaCRec - estimate of ThetaCerenkov
217 //              status    - TRUE if a good solution is found, FALSE if not
218
219   Double_t *phiphot = new Double_t[fNClu];  
220   Double_t *dist    = new Double_t[fNClu];  
221   Int_t    *indphi  = new    Int_t[fNClu];  
222
223   Bool_t status;
224     
225 // Sort in phi angle...  
226   for(Int_t i=0;i<fNClu;i++) {
227     phiphot[i] = (TMath::ATan2(fMipY-fYClu[i],fMipX-fXClu[i])+TMath::Pi())*TMath::RadToDeg();
228     dist[i]=TMath::Sqrt((fMipX-fXClu[i])*(fMipX-fXClu[i])+(fMipY-fYClu[i])*(fMipY-fYClu[i]));
229   }
230   
231   TMath::Sort(fNClu,phiphot,indphi,kFALSE);
232   
233 // Purify with a truncated mean;
234   Int_t np=0;
235   Double_t dMean  = 0;
236   Double_t dMean2 = 0;
237   for(Int_t i=0;i<fNClu;i++) {
238     if(!fClCk[indphi[i]]) continue;                                                  // Check if a good photon candidate or not
239     dMean +=dist[indphi[i]];
240     dMean2+=dist[indphi[i]]*dist[indphi[i]];
241     np++;
242   }
243   
244   dMean  /=(Double_t)np;
245   dMean2 /=(Double_t)np;
246   Double_t rms = TMath::Sqrt(dMean2 - dMean*dMean);
247   
248   for(Int_t i=0;i<fNClu;i++) {
249     if(!fClCk[indphi[i]]) continue;                                                  // Check if a good photon candidate or not
250     if(TMath::Abs(dMean-dist[indphi[i]]) > 2*rms) {
251       fClCk[indphi[i]] = kFALSE;
252       continue;
253     }
254   }
255
256 //
257 //  purify vectors for good photon candidates
258 //
259   Int_t npeff=0;
260   Double_t *phiphotP = new Double_t[fNClu+1];  
261   Double_t *distP    = new Double_t[fNClu+1];  
262   for(Int_t i=0;i<fNClu;i++) {
263     if(!fClCk[indphi[i]]) continue;                                                  // Check if a good photon candidate or not
264     phiphotP[npeff] = phiphot[indphi[i]];
265     distP[npeff]    = dist[indphi[i]];
266     npeff++;
267   }
268   
269   delete [] phiphot;
270   delete [] dist;
271   delete [] indphi;
272
273   if(npeff<3) {
274     delete [] phiphotP;
275     delete [] distP;
276     return kFALSE;
277   }
278
279 //  for(Int_t i=0;i<npeff;i++) {Printf(" n. %d phiphot %f dist %f",i,phiphotP[i],distP[i]);}
280   
281   Double_t xA,xB;
282   if(ShapeModel(npeff,phiphotP,distP,xA,xB,phiTrkRec)) {
283     
284 //    Printf("FindShape: phi start %f",phiTrkRec*TMath::RadToDeg());
285
286     Int_t bin = fgDatabase->FindBin(xA,xB);
287     Int_t compact = (Int_t)fgDatabase->GetBinContent(bin);
288     thetaCRec =  (Double_t)(compact%1000);
289     thTrkRec  =  (Double_t)(compact/1000);
290
291     thTrkRec *= TMath::DegToRad(); 
292     thetaCRec *= TMath::DegToRad();
293
294 //    Printf("FindShape: xA %f xB %f compact %d thTrk %f thC %f",xA,xB,compact,thTrkRec*TMath::RadToDeg(),thetaCRec*TMath::RadToDeg());
295         
296     status = kTRUE;
297     
298   } else {
299     
300     status = kFALSE;
301     
302   }
303
304   delete [] phiphotP;
305   delete [] distP;
306   
307   return status;
308 }
309 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
310 Bool_t AliHMPIDReconHTA::ShapeModel(Int_t np,Double_t *phiphot,Double_t *dist,Double_t &xA,Double_t &xB,Double_t &phiStart)
311 {
312 // Find a Spline curve to define dist. vs. phi angle
313 // in order to estimate the phi of the track
314 // Arguments:   np     - # points corresponding to # photon candidates
315 //             dist    - distance of each photon from MIP
316 //             phiphot - phi of the photon in the DRS
317 //   Returns:  xA      - min. distance from MIP
318 //             xB      - dist. from mip perpedicular to the major axis 
319 //             phiStart- estimate of the track phi
320
321   TGraph *phigr = new TGraph(np,phiphot,dist);
322   TSpline3 *sphi = new TSpline3("sphi",phigr);
323   if(!sphi) {Printf("Spline not created!Bye.");return kFALSE;}
324
325   Int_t locMin = TMath::LocMin(np,dist);
326   Int_t locMax = TMath::LocMax(np,dist);
327   
328   Double_t minX = phiphot[locMin];
329 //  Double_t minY =    dist[locMin];
330   Double_t maxX = phiphot[locMax];
331 //  Double_t maxY =    dist[locMax];
332   
333   Int_t ip[3] = {-1,0,1};
334   if(locMin==0   ) {ip[0]= 0;ip[1]= 1;ip[2]= 2;}
335   if(locMin==np-1) {ip[0]=-2;ip[1]=-1;ip[2]= 0;}
336   
337   Double_t minXf = VertParab(phiphot[locMin+ip[0]],dist[locMin+ip[0]], 
338                              phiphot[locMin+ip[1]],dist[locMin+ip[1]],
339                              phiphot[locMin+ip[2]],dist[locMin+ip[2]]);
340   if(minXf< phiphot[locMin+ip[0]] || minXf > phiphot[locMin+ip[2]]) minXf = minX;
341   
342   ip[0]=-1;ip[1]=0;ip[2]=1;
343   if(locMax==0   ) {ip[0]= 0;ip[1]= 1;ip[2]= 2;}
344   if(locMax==np-1) {ip[0]=-2;ip[1]=-1;ip[2]= 0;}
345   
346   Double_t maxXf = VertParab(phiphot[locMax+ip[0]],dist[locMax+ip[0]], 
347                              phiphot[locMax+ip[1]],dist[locMax+ip[1]],
348                              phiphot[locMax+ip[2]],dist[locMax+ip[2]]);
349   if(maxXf< phiphot[locMax+ip[0]] || maxXf > phiphot[locMax+ip[2]]) maxXf = maxX;
350   
351 //  Printf(" phi at mindist %f and found %f",minX,minXf);
352 //  Printf(" phi at maxdist %f and found %f",maxX,maxXf);
353 //  
354   if(TMath::Abs(maxXf-minXf)>30) {
355     xA = sphi->Eval(minXf);
356     if(xA < 0) xA = dist[sphi->FindX(xA)];
357     xB = sphi->Eval(minXf-90);
358     if(xB < 0) xB = dist[sphi->FindX(xB)];
359     phiStart = minXf-180;  //open ring or acceptance effect...so believe to min phi angle!
360   } else {
361     phiStart = 0.5*(maxXf-180+minXf);
362     xA = sphi->Eval(phiStart+180);
363     if(xA < 0) xA = dist[sphi->FindX(xA)];
364     xB = sphi->Eval(phiStart+90);
365     if(xB < 0) xB = dist[sphi->FindX(xB)];
366   }
367   //
368 //  Printf("ShapeModel: phiStart %f xA %f xB %f",phiStart,xA,xB);
369   
370   phiStart*=TMath::DegToRad();
371   return kTRUE;
372 }
373 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
374 Double_t AliHMPIDReconHTA::VertParab(Double_t x1,Double_t y1,Double_t x2, Double_t y2, Double_t x3, Double_t y3)const
375 {
376 // It uses parabola from 3 points to evaluate the x-coord of the parab 
377 // Arguments:    xi,yi - points
378 //   Returns:    x-coord of the vertex 
379   
380   Double_t a = ((x1-x3)*(y1-y2)-(x1-x2)*(y1-y3))/((x1*x1-x2*x2)*(x1-x3)-(x1*x1-x3*x3)*(x1-x2));
381   Double_t b = (y1-y2 - a*(x1*x1-x2*x2))/(x1-x2);
382 //  Double_t c = y1 - a*x1*x1-b*x1;
383   return -0.5*b/a;
384 }
385 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
386 Bool_t AliHMPIDReconHTA::FitFree(Double_t thTrkRec,Double_t phiTrkRec)
387 {
388 // Fit performed by minimizing RMS/sqrt(n) of the
389 // photons reconstructed. First phi is fixed and theta
390 // is fouond, then (th,ph) of the track
391 // as free parameters
392 // Arguments:    PhiRec phi of the track
393 //   Returns:    none
394   
395   TMinuit *pMinuit = new TMinuit(2);
396   pMinuit->mncler();
397   gMinuit->SetObjectFit((TObject*)this);  gMinuit->SetFCN(AliHMPIDReconHTA::FunMinPhot);  //set fit function
398   Double_t aArg=-1,parStep,parLow,parHigh;     Int_t iErrFlg;                 //tmp vars for TMinuit
399   Double_t d1,d2,d3;
400   TString sName;
401   Double_t th,ph;
402   
403   pMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg);                                          //suspend all printout from TMinuit
404   pMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg);
405
406   if(thTrkRec==0) thTrkRec = 3.*TMath::DegToRad();    // not to start from the edge...
407   
408   pMinuit->mnparm(0," thTrk  ",thTrkRec ,parStep=0.01,parLow=0,parHigh=TMath::PiOver4(),iErrFlg);
409   pMinuit->mnparm(1," phiTrk ",phiTrkRec,parStep=0.01,parLow=0,parHigh=TMath::TwoPi(),iErrFlg);
410   
411   pMinuit->FixParameter(1);
412   pMinuit->mnexcm("SIMPLEX" ,&aArg,0,iErrFlg);   
413   pMinuit->mnexcm("MIGRAD"  ,&aArg,0,iErrFlg);
414   pMinuit->Release(1);  
415   pMinuit->mnexcm("MIGRAD"  ,&aArg,0,iErrFlg);
416   
417   pMinuit->mnpout(0,sName,th,d1,d2,d3,iErrFlg);
418   pMinuit->mnpout(1,sName,ph,d1,d2,d3,iErrFlg);   
419   
420   Double_t f,par[2];
421   Double_t *grad=0x0;
422   par[0] = th;par[1] = ph;
423   pMinuit->Eval(2,grad,f,par,3);
424
425 //  Printf("FitFree: theta %f phi %f",th,ph);
426   
427   SetTrkFit(th,ph);
428   return kTRUE;
429 }
430 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
431 void AliHMPIDReconHTA::FunMinPhot(Int_t &/* */,Double_t* /* */,Double_t &f,Double_t *par,Int_t iflag)
432 {
433 // Minimization function to find best track and thetaC parameters
434 // Arguments:    f = function value to minimize
435 //             par = list of parameter to find
436 //           iflag = flag status. See Minuit instructions
437 //   Returns:    none
438 //
439 // Note: it is necessary to call an instance of AlihMPIDParam. Not possible to use fParam
440 // because of the static instantiation of the function in Minuit
441   
442   AliHMPIDParam *pParam=AliHMPIDParam::Instance();
443   AliHMPIDReconHTA *pRecHTA=(AliHMPIDReconHTA*)gMinuit->GetObjectFit();
444   AliHMPIDRecon pRec;
445   Double_t sizeCh = 0.5*pParam->RadThick()+pParam->WinThick()+pParam->GapThick();
446   Double_t thTrk = par[0]; 
447   Double_t phTrk = par[1];
448   Double_t xrad = pRecHTA->MipX() - sizeCh*TMath::Tan(thTrk)*TMath::Cos(phTrk);
449   Double_t yrad = pRecHTA->MipY() - sizeCh*TMath::Tan(thTrk)*TMath::Sin(phTrk);
450   pRecHTA->SetRadXY(xrad,yrad);
451   pRec.SetTrack(xrad,yrad,thTrk,phTrk);
452
453   Double_t meanCkov =0;
454   Double_t meanCkov2=0;
455   Double_t thetaCer,phiCer;
456   Int_t nClAcc = 0;
457   Int_t nClTot=pRecHTA->NClu();
458     
459   for(Int_t i=0;i<nClTot;i++) {
460     if(!(pRecHTA->ClCk(i))) continue;
461     pRec.FindPhotCkov(pRecHTA->XClu(i),pRecHTA->YClu(i),thetaCer,phiCer);  
462     meanCkov  += thetaCer;
463     meanCkov2 += thetaCer*thetaCer;
464     nClAcc++;
465   }
466   if(nClAcc==0) {f=999;return;}
467   meanCkov  /=(Double_t)nClAcc;
468   meanCkov2 /=(Double_t)nClAcc;
469   Double_t rms = TMath::Sqrt(TMath::Abs(meanCkov2 - meanCkov*meanCkov));
470   f = rms/TMath::Sqrt((Double_t)nClAcc);
471   
472   if(iflag==3) {
473 /*
474     Printf("FunMinPhot before: photons candidates %i used %i",nClTot,nClAcc);
475     nClAcc = 0;
476     Double_t meanCkov1=0;
477     Double_t meanCkov2=0;
478     for(Int_t i=0;i<nClTot;i++) {
479       if(!(pRec->ClCk(i))) continue;
480       pRec->FindPhotCkov(pRec->XClu(i),pRec->YClu(i),thetaCer,phiCer);  
481       if(TMath::Abs(thetaCer-meanCkov)<2*rms) {
482         meanCkov1 += thetaCer;
483         meanCkov2 += thetaCer*thetaCer;
484         nClAcc++;
485       } else pRec->SetClCk(i,kFALSE);
486     }
487     meanCkov1/=nClAcc;
488     Double_t rms2 = (meanCkov2 - meanCkov*meanCkov*nClAcc)/nClAcc;
489     Printf("FunMinPhot after: photons candidates %i used %i thetaCer %f",nClTot,nClAcc,meanCkov1);
490     pRec->SetCkovFit(meanCkov1);
491     pRec->SetCkovSig2(rms2);
492     pRec->SetNClu(nClAcc);
493 */
494 //    Printf("FunMinPhot: reconstructed theta Cerenkov %f with %d photons",meanCkov,nClAcc);
495     pRecHTA->SetCkovFit(meanCkov);
496     pRecHTA->SetCkovSig2(rms*rms);
497     pRecHTA->SetNClu(nClAcc);
498   }
499 }//FunMinPhot()
500 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
501 void AliHMPIDReconHTA::InitDatabase()
502 {
503 // Construction a database of ring shapes on fly
504 //   Arguments: none
505 //   Returns  : none
506 //  N.B. fgDatabase points to a TH2I with x-min dist from MIP
507 //                                        y-dist from the ring of the MIP perpendicular to major axis
508 //        The content is the packed info of track theta and thetaC in degrees
509 //                        thetaC+1000*thTrk
510 //
511   TStopwatch timer;
512   timer.Start();
513   
514   AliInfo(Form("database HTA is being built.Please, wait..."));
515 //  
516   Double_t x[3]={0,0,0},y[3];
517
518   AliHMPIDRecon rec;
519
520   if(!fParam) fParam=AliHMPIDParam::Instance();
521   Double_t thetaMax = TMath::ACos(1./fParam->MeanIdxRad());
522   Double_t thTrkMax = 1./TMath::ASin(fParam->MeanIdxRad());    
523
524   Int_t nstepx = 1000;
525   Int_t nstepy = 1000;
526
527   TH2I *deconv = new TH2I("deconv","database;d1;d2;thC+1000*thTrk",500,0,50,150,0,15);
528   //
529   Double_t xrad = 0;
530   Double_t yrad = 0;
531   Double_t phTrk = 0;
532
533   for(Int_t i=0;i<nstepx;i++) {     //loop on thetaC
534     for(Int_t j=0;j<nstepy;j++) {   //loop on theta particle
535       Double_t thetaC = thetaMax/nstepx*((Double_t)i+0.5);
536       Double_t thTrk  = thTrkMax/nstepy*((Double_t)j+0.5);
537       //
538       //mip position
539       //
540       Double_t sizeCh = fParam->RadThick()+fParam->WinThick()+fParam->GapThick();
541       Double_t xmip = xrad + sizeCh*TMath::Tan(thTrk)*TMath::Cos(phTrk);
542       Double_t ymip = yrad + sizeCh*TMath::Tan(thTrk)*TMath::Sin(phTrk);
543
544       Double_t dist1,dist2;
545       //
546       //first point at phi=0
547       //
548       rec.SetTrack(xrad,yrad,thTrk,phTrk);
549       TVector2 pos;
550       pos=rec.TracePhot(thetaC,0);
551
552       if(pos.X()==-999) {
553         dist1 = 0;            //open ring...anly the distance btw mip and point at 180 will be considered
554       } else {
555         x[0] = pos.X(); y[0] = pos.Y();
556         dist1   = TMath::Sqrt((x[0]-xmip)*(x[0]-xmip)+(y[0]-ymip)*(y[0]-ymip));
557       }
558       //
559       //second point at phi=180
560       //
561       rec.SetTrack(xrad,yrad,thTrk,phTrk);
562       pos=rec.TracePhot(thetaC,TMath::Pi());
563
564       if(pos.X()==-999) {Printf("it should not happens!Bye");return;}
565       x[1] = pos.X(); y[1] = pos.Y();
566       if((x[1]-xmip)*(x[0]-xmip)>0) continue; // to avoid circles out mips (for very low ThetaC)
567       dist2   = TMath::Sqrt((x[1]-xmip)*(x[1]-xmip)+(y[1]-ymip)*(y[1]-ymip));
568
569 //      Double_t distA = dist1+dist2;
570       Double_t distA = dist2;     // only the minimum: problem of acceptance
571       //
572       //second point at phi=90
573       //
574       rec.SetTrack(xrad,yrad,thTrk,phTrk);
575       pos=rec.TracePhot(thetaC,TMath::PiOver2());
576
577       if(pos.X()==-999) continue;
578       x[2] = pos.X(); y[2] = pos.Y();
579       Double_t distB   = TMath::Sqrt((x[2]-xmip)*(x[2]-xmip)+(y[2]-ymip)*(y[2]-ymip));
580 // compact the infos...      
581       Int_t compact = (Int_t)(thetaC*TMath::RadToDeg())+1000*(Int_t)(thTrk*TMath::RadToDeg());
582       Int_t bin = deconv->FindBin(distA,distB);
583       if(deconv->GetBinContent(bin)==0) deconv->Fill(distA,distB,compact);
584     }
585   }
586
587   FillZeroChan(deconv);
588   fgDatabase = deconv;
589
590   timer.Stop();
591   Double_t nSecs = timer.CpuTime();  
592   AliInfo(Form("database HTA successfully open in %3.1f sec.(CPU). Reconstruction is started.",nSecs));
593 }
594 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
595 void AliHMPIDReconHTA::FillZeroChan(TH2I *deconv)const
596 {
597   //If fills eventually channel without entries
598   //inthe histo "database" jyst interpolating the neighboring cells
599   // Arguments: histogram pointer of the database
600   //   Returns: none
601   //
602   Int_t nbinx = deconv->GetNbinsX();
603   Int_t nbiny = deconv->GetNbinsY();
604   for(Int_t i = 0;i<nbinx;i++) {
605     for(Int_t j = 0;j<nbiny;j++) {
606       if(deconv->GetBinContent(i,j) == 0) {
607         Int_t nXmin = i-1; Int_t nXmax=i+1;
608         Int_t nYmin = j-1; Int_t nYmax=j+1;
609         Int_t nc = 0;
610         Double_t meanC =0;
611         Double_t meanTrk =0;
612         for(Int_t ix=nXmin;ix<=nXmax;ix++) {
613           if(ix<0||ix>nbinx) continue;
614           for(Int_t iy=nYmin;iy<=nYmax;iy++) {
615             if(iy<0||iy>nbiny) continue;
616             meanC  +=  (Int_t)deconv->GetBinContent(ix,iy)%1000;
617             meanTrk+=  (Int_t)deconv->GetBinContent(ix,iy)/1000;
618             nc++;
619           }
620           meanC/=nc; meanTrk/=nc;
621           Int_t compact = (Int_t)meanC+1000*(Int_t)meanTrk;
622           if(compact>0)deconv->SetCellContent(i,j,compact);
623         }
624       }
625     }
626   }
627 }
628 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++