]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDTracker.cxx
TH1 fitting protections improved
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDTracker.cxx
CommitLineData
3c6274c1 1#include "AliHMPIDTracker.h" //class header
e56b695f 2#include "AliHMPIDtrack.h" //class header
3c6274c1 3#include "AliHMPIDCluster.h" //GetTrackPoint(),PropagateBack()
4#include "AliHMPIDParam.h" //GetTrackPoint(),PropagateBack()
e56b695f 5#include "AliHMPIDPid.h" //Recon(),reconHTA()
59d9d4b3 6#include "AliHMPIDRecon.h" //Recon()
c458aab6 7#include "AliHMPIDRecoParamV1.h" //Recon()
b8541917 8#include "AliHMPIDReconstructor.h"//Recon()
57fcdc2b 9#include "AliHMPIDReconHTA.h" //ReconHTA()
10#include <AliESDEvent.h> //PropagateBack(),Recon()
c38d443f 11#include <AliESDtrack.h> //Intersect()
12#include <AliTracker.h>
d3da6dc4 13#include <AliRun.h> //GetTrackPoint(),PropagateBack()
14#include <AliTrackPointArray.h> //GetTrackPoint()
15#include <AliAlignObj.h> //GetTrackPoint()
59d9d4b3 16#include <AliCDBManager.h> //PropageteBack()
17#include <AliCDBEntry.h> //PropageteBack()
423554a3 18//.
19// HMPID base class fo tracking
20//.
21//.
22//.
d3da6dc4 23ClassImp(AliHMPIDTracker)
94b1fbfa 24//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
c61a7285 25AliHMPIDTracker::AliHMPIDTracker():
26 AliTracker(),
27 fClu(new TObjArray(AliHMPIDParam::kMaxCh+1))
94b1fbfa 28{
59d9d4b3 29// ctor. Create TObjArray of TClonesArray of AliHMPIDCluster
30//
31//
c61a7285 32 fClu->SetOwner(kTRUE);
ae5a42aa 33 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
59d9d4b3 34}//ctor
d3da6dc4 35//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
36Bool_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()
59d9d4b3 44 Int_t iCham=idx/1000000; Int_t iClu=idx%1000000;
2247c219 45 iClu = iClu%1000; //GetHMPIDcluIdx -> 1e+6*ch + 1e+3*clusize + cluIdx;
46 point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham));//layer and chamber number
94b1fbfa 47 TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham];
48 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
2247c219 49 Float_t xyz[3];
50 pClu->GetGlobalXYZ(xyz);
51 Float_t cov[6];
52 pClu->GetGlobalCov(cov);
53 point.SetXYZ(xyz,cov);
d3da6dc4 54 return kTRUE;
59d9d4b3 55}//GetTrackPoint()
56//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
39cd22e6 57Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
59d9d4b3 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
b4866751 62 AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching
ae5a42aa 63 for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){ //chambers loop
b4866751 64 Int_t chInt = IntTrkCha(i,hmpTrk,xPc,yPc,xRa,yRa,theta,phi);
65 if(chInt>=0) {delete hmpTrk;hmpTrk=0x0;return chInt;}
59d9d4b3 66 } //chambers loop
b4866751 67 delete hmpTrk; hmpTrk=0x0;
39cd22e6 68 return -1; //no intersection with HMPID chambers
59d9d4b3 69}//IntTrkCha()
d3da6dc4 70//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
496c71b0 71Int_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];
b4866751 81 pParam->Norm(ch,n2);
496c71b0 82 pParam->Point(ch,p2,AliHMPIDParam::kPc); //point & norm for entrance to PC plane
97a4d538 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;
496c71b0 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//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
d3da6dc4 92Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
93{
59d9d4b3 94// Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
d3da6dc4 95// Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
96// Returns: error code (currently ignored in AliReconstruction::RunTraking())
ae5a42aa 97 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
94b1fbfa 98 pCluTree->GetEntry(0);
94b1fbfa 99 return 0;
59d9d4b3 100}//LoadClusters()
d3da6dc4 101//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
af885e0f 102Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
abb5f786 103{
59d9d4b3 104// Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
abb5f786 105// Agruments: pEsd - pointer to ESD
106// Returns: error code
f455af6e 107 AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
49881df7 108 AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
0a7c7084 109 if(!pNmeanEnt) AliError("No Nmean C6F14 ");
110 if(!pQthreEnt) AliError("No Qthre");
94b1fbfa 111
afe12692 112 return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());
611e810d 113}//PropagateBack()
abb5f786 114//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
afe12692 115Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
d3da6dc4 116{
59d9d4b3 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)
d3da6dc4 119// Returns: error code, 0 if no errors
496c71b0 120
121 AliHMPIDRecon recon; //instance of reconstruction class, nothing important in ctor
b8541917 122 AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam
39cd22e6 123 Float_t xPc,yPc,xRa,yRa,theta,phi;
43e51ec6 124 Double_t cluLORS[2]={0};
125// Double_t cluMARS[3]={0},trkMARS[3]={0};
496c71b0 126// Double_t bestcluMARS[3]={0,0,0};
43e51ec6 127// Double_t radClu,radInitTrk;
496c71b0 128 Int_t nMipClusTot=0;
43e51ec6 129// Double_t d3d=0;
84770d7b 130 Double_t qthre = 0; Double_t nmean=0; Int_t hvsec=0;
496c71b0 131 Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
b8541917 132
496c71b0 133
134 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event
c38d443f 135
43e51ec6 136// Double_t bestChi2=99999;chi2=99999; //init. track matching params
b8541917 137 Double_t dmin=999999,bz=0,distCut=1,distParams[5]={1};
43e51ec6 138
139 Bool_t isOkDcut=kFALSE;
140 Bool_t isOkQcut=kFALSE;
141 Bool_t isMatched=kFALSE;
142
496c71b0 143 AliHMPIDCluster *bestHmpCluster=0x0; //the best matching cluster
c38d443f 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
496c71b0 150 AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching
151 bz=AliTracker::GetBz();
43e51ec6 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
496c71b0 157
af291e40 158 Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi); //find the intersected chamber for this track
43e51ec6 159 if(ipCh<0) continue; //no intersection at all, go after next track
160
97a4d538 161 pTrk->SetHMPIDtrk(xPc,yPc,theta,phi); //store initial infos
43e51ec6 162 pTrk->SetHMPIDcluIdx(ipCh,9999); //set chamber, index of cluster + cluster size
496c71b0 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
43e51ec6 170 if(nMipClusTot==0) continue;
171
172 Int_t index=-1; //index of the "best" matching cluster
173
496c71b0 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) { // just for backward compatibility
179 qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(pEsd->GetTimeStamp()); //
180 } else { // in the past just 1 qthre
181 hvsec = pParam->InHVSector(pClu->Y()); // per chamber
182 if(hvsec>=0)
84770d7b 183 qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(pEsd->GetTimeStamp()); //
496c71b0 184 } //
185//
186 if(pClu->Q()<qthre) continue; //charge compartible with MIP clusters
187 isOkQcut = kTRUE;
43e51ec6 188//
496c71b0 189 cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y(); //get the LORS coordinates of the cluster
43e51ec6 190 Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
191
192 if(dist<dmin) {
193 dmin = dist;
194 index=iClu;
195 bestHmpCluster=pClu;
196 }
197 }
b8541917 198
496c71b0 199 if(!isOkQcut) {
496c71b0 200 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
201 continue;
202 }
c38d443f 203
204 Double_t radius = (pParam->Lors2Mars(ipCh,pParam->SizeAllX()/2,pParam->SizeAllY()/2)).Mag();
205
87a6d94b 206 if(!AliTracker::PropagateTrackToBxByBz(hmpTrk,radius,pTrk->GetMass(),1,kFALSE)) continue;
43e51ec6 207
c38d443f 208 if(!hmpTrk->PropagateTo(bestHmpCluster)) continue;
209
43e51ec6 210 Int_t cluSiz = bestHmpCluster->Size();
211 pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0); //store mip info in any case
212 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set chamber, index of cluster + cluster size
b8541917 213
214 if(AliHMPIDReconstructor::GetRecoParam()) //retrieve distance cut
215 {
b8541917 216 if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE) //distance cut is fixed number
217 {
218 distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist();
219 }
220 else
221 {
13fb7604 222 for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar); //prevision: distance cut is function of momentum
b8541917 223 distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise
224 }
225 }
226 else {distCut=pParam->DistCut();}
227
228 if(dmin < distCut) {
496c71b0 229 isOkDcut = kTRUE;
230 }
231
b8541917 232
233
496c71b0 234 if(!isOkDcut) {
496c71b0 235 pTrk->SetHMPIDsignal(pParam->kMipDistCut); //closest cluster with enough charge is still too far from intersection
236 }
237
238 if(isOkQcut*isOkDcut) isMatched = kTRUE; // MIP-Track matched !!
239
240 if(!isMatched) continue; // If matched continue...
43e51ec6 241
242 Bool_t isOk = hmpTrk->Update(bestHmpCluster,0.1,0);
243 if(!isOk) continue;
c38d443f 244 pTrk->SetOuterHmpParam(hmpTrk,AliESDtrack::kHMPIDout);
245
246 FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
b8541917 247
496c71b0 248
496c71b0 249 //evaluate nMean
250 if(pNmean->GetEntries()==21) { //for backward compatibility
84770d7b 251 nmean=((TF1*)pNmean->At(3*ipCh))->Eval(pEsd->GetTimeStamp()); //C6F14 Nmean for this chamber
496c71b0 252 } else {
253 Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved
0a7c7084 254 if(iRad < 0) {
255 nmean = -1;
256 } else {
84770d7b 257 Double_t tLow = ((TF1*)pNmean->At(6*ipCh+2*iRad ))->Eval(pEsd->GetTimeStamp()); //C6F14 low temp for this chamber
258 Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(pEsd->GetTimeStamp()); //C6F14 high temp for this chamber
496c71b0 259 Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y
260 nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp
0a7c7084 261 }
496c71b0 262 if(nmean < 0){ //track didn' t pass through the radiator
263 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag
264 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set index of cluster
265 continue;
266 }
267 }
268 //
a591e55f 269 recon.SetImpPC(xPc,yPc); //store track impact to PC
97a4d538 270 recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa); //search for Cerenkov angle of this track
235a21c2 271
272 if(pTrk->GetHMPIDsignal()<0) continue;
273
e56b695f 274 AliHMPIDPid pID;
275 Double_t prob[5];
276 pID.FindPid(pTrk,5,prob);
277 pTrk->SetHMPIDpid(prob);
278// 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);
279
b4866751 280 delete hmpTrk;hmpTrk=0x0;
496c71b0 281 }//iTrk
282
d3da6dc4 283 return 0; // error code: 0=no error;
611e810d 284}//Recon()
285//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
9785d5fb 286Int_t AliHMPIDTracker::ReconHiddenTrk(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
611e810d 287{
288// Static method to reconstruct Theta Ckov for all valid tracks of a given event.
afe12692 289// 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)
611e810d 290// Returns: error code, 0 if no errors
9785d5fb 291
292 AliHMPIDReconHTA reconHTA; //instance of reconstruction class, nothing important in ctor
293
294 AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam
295
296 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event
297
298 AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //here it is simulated or just empty track
299 Int_t ipCh;
300 ipCh = pTrk->GetHMPIDcluIdx();ipCh/=1000000;
301 if(ipCh<0) continue;
302
303 TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters
304 Int_t nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber
305
306 Double_t qMip=-1;
307 Int_t chMip=-1;
ee2f3539 308 Double_t xMip = 0;
309 Double_t yMip = 0;
310 Int_t indMip = 0;
311 Int_t cluMipSiz = 0;
9785d5fb 312
313 for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop
314
315 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster
316 Double_t qClus = pClu->Q();
317 if(qClus>qMip) {
318 qMip = qClus;
319 chMip = pClu->Ch();
320 xMip = pClu->X();
321 yMip = pClu->Y();
322 indMip = iClu;
323 cluMipSiz = pClu->Size();
324 }
325 }//clus loop
326
327 if(chMip<0) return 1;
328
329 Int_t hvsec;
ee2f3539 330 Double_t qthre=0;
9785d5fb 331// evaluate qThre
332 if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) { // just for backward compatibility
333 qthre=((TF1*)pQthre->At(chMip))->Eval(pEsd->GetTimeStamp()); //
334 } else { // in the past just 1 qthre
335 hvsec = pParam->InHVSector(yMip); // per chamber
336 if(hvsec>=0) qthre=((TF1*)pQthre->At(6*chMip+hvsec))->Eval(pEsd->GetTimeStamp()); //
337 }
338//
339 if(qMip<qthre) {
340 pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case
341 pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);
342 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
343 return 1; //charge compatible with MIP clusters
344 }
345
346 pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case
347 pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);
348
349 Double_t yRa = yMip; //just an approx...
350 Double_t nmean;
351 //evaluate nMean
352 if(pNmean->GetEntries()==21) { //for backward compatibility
353 nmean=((TF1*)pNmean->At(3*chMip))->Eval(pEsd->GetTimeStamp()); //C6F14 Nmean for this chamber
354 } else {
355 Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved
356 if(iRad < 0) {
357 nmean = -1;
358 } else {
359 Double_t tLow = ((TF1*)pNmean->At(6*chMip+2*iRad ))->Eval(pEsd->GetTimeStamp()); //C6F14 low temp for this chamber
360 Double_t tHigh = ((TF1*)pNmean->At(6*chMip+2*iRad+1))->Eval(pEsd->GetTimeStamp()); //C6F14 high temp for this chamber
361 Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y
362 nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp
363 }
364 if(nmean < 0){ //track didn' t pass through the radiator
365 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag
366 return 1;
367 }
368 }
369 //
ee2f3539 370 if(!reconHTA.CkovHiddenTrk(pTrk,(TClonesArray *)pClus->At(ipCh),indMip,nmean)) { //search for track parameters and Cerenkov angle of this track
371 AliHMPIDPid pID;
372 Double_t prob[5];
373 pID.FindPid(pTrk,5,prob);
374 pTrk->SetHMPIDpid(prob);
375 }
e56b695f 376// 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);
9785d5fb 377 }//iTrk
378
379 return 0; // error code: 0=no error;
380
611e810d 381}//Recon()
d3da6dc4 382//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
c1af14f7 383void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
384
385 // Publishes all pointers to clusters known to the tracker into the
386 // passed object array.
387 // The ownership is not transfered - the caller is not expected to delete
388 // the clusters
389
390 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
391 TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
392 for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
393 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);
394 array->AddLast(pClu);
395 }//cluster loop in iCh
396 pCluArr->Delete();
397 }//Ch loop
398
399 return;
400}
401//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++