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