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