]>
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() |
b4b2d95f | 7 | #include "AliHMPIDRecoParamV1.h" //Recon() |
b8541917 | 8 | #include "AliHMPIDReconstructor.h"//Recon() |
57fcdc2b | 9 | #include "AliHMPIDReconHTA.h" //ReconHTA() |
3e630437 | 10 | #include <AliLog.h> //Recon() |
57fcdc2b | 11 | #include <AliESDEvent.h> //PropagateBack(),Recon() |
c38d443f | 12 | #include <AliESDtrack.h> //Intersect() |
13 | #include <AliTracker.h> | |
d3da6dc4 | 14 | #include <AliRun.h> //GetTrackPoint(),PropagateBack() |
15 | #include <AliTrackPointArray.h> //GetTrackPoint() | |
16 | #include <AliAlignObj.h> //GetTrackPoint() | |
59d9d4b3 | 17 | #include <AliCDBManager.h> //PropageteBack() |
18 | #include <AliCDBEntry.h> //PropageteBack() | |
79439569 | 19 | #include "TTreeStream.h" // debug streamer |
20 | // | |
423554a3 | 21 | // HMPID base class fo tracking |
22 | //. | |
23 | //. | |
24 | //. | |
d3da6dc4 | 25 | ClassImp(AliHMPIDTracker) |
94b1fbfa | 26 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
c61a7285 | 27 | AliHMPIDTracker::AliHMPIDTracker(): |
28 | AliTracker(), | |
79439569 | 29 | fClu(new TObjArray(AliHMPIDParam::kMaxCh+1)), |
30 | fDebugStreamer(0) | |
94b1fbfa | 31 | { |
59d9d4b3 | 32 | // ctor. Create TObjArray of TClonesArray of AliHMPIDCluster |
33 | // | |
34 | // | |
c61a7285 | 35 | fClu->SetOwner(kTRUE); |
ae5a42aa | 36 | for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i); |
79439569 | 37 | fDebugStreamer = new TTreeSRedirector("HMPIDdebug.root"); |
38 | ||
59d9d4b3 | 39 | }//ctor |
79439569 | 40 | |
41 | ||
42 | AliHMPIDTracker::~AliHMPIDTracker(){ | |
43 | // | |
44 | // destructor | |
45 | // | |
46 | delete fClu; | |
47 | if (fDebugStreamer) delete fDebugStreamer; | |
48 | } | |
49 | ||
d3da6dc4 | 50 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
51 | Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const | |
52 | { | |
53 | // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track. | |
54 | // MIP cluster is reffered by index which is stored in AliESDtrack ??????? | |
55 | // Arguments: idx- cluster index which is stored by HMPID in AliESDtrack | |
56 | // point- reference to the object where to store the point | |
57 | // Returns: status of operation if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. | |
58 | if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack() | |
59d9d4b3 | 59 | Int_t iCham=idx/1000000; Int_t iClu=idx%1000000; |
2247c219 | 60 | iClu = iClu%1000; //GetHMPIDcluIdx -> 1e+6*ch + 1e+3*clusize + cluIdx; |
61 | point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham));//layer and chamber number | |
94b1fbfa | 62 | TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham]; |
53b3604e | 63 | if(!pArr) return kFALSE; |
94b1fbfa | 64 | AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster |
53b3604e | 65 | if(!pClu) return kFALSE; |
2247c219 | 66 | Float_t xyz[3]; |
67 | pClu->GetGlobalXYZ(xyz); | |
68 | Float_t cov[6]; | |
69 | pClu->GetGlobalCov(cov); | |
70 | point.SetXYZ(xyz,cov); | |
d3da6dc4 | 71 | return kTRUE; |
59d9d4b3 | 72 | }//GetTrackPoint() |
73 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
39cd22e6 | 74 | 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 | 75 | { |
76 | // Static method to find intersection in between given track and HMPID chambers | |
77 | // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm] | |
78 | // Returns: intersected chamber ID or -1 | |
b4866751 | 79 | AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching |
ae5a42aa | 80 | for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){ //chambers loop |
b4866751 | 81 | Int_t chInt = IntTrkCha(i,hmpTrk,xPc,yPc,xRa,yRa,theta,phi); |
82 | if(chInt>=0) {delete hmpTrk;hmpTrk=0x0;return chInt;} | |
59d9d4b3 | 83 | } //chambers loop |
b4866751 | 84 | delete hmpTrk; hmpTrk=0x0; |
39cd22e6 | 85 | return -1; //no intersection with HMPID chambers |
59d9d4b3 | 86 | }//IntTrkCha() |
d3da6dc4 | 87 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
496c71b0 | 88 | 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) |
89 | { | |
90 | // Static method to find intersection in between given track and HMPID chambers | |
91 | // Arguments: pTrk- HMPID track; xPc,yPc- track intersection with PC in LORS [cm] | |
92 | // Returns: intersected chamber ID or -1 | |
93 | AliHMPIDParam *pParam=AliHMPIDParam::Instance(); | |
94 | Double_t p1[3],n1[3]; | |
95 | pParam->Norm(ch,n1); | |
96 | pParam->Point(ch,p1,AliHMPIDParam::kRad); //point & norm for middle of radiator plane | |
97 | Double_t p2[3],n2[3]; | |
b4866751 | 98 | pParam->Norm(ch,n2); |
496c71b0 | 99 | pParam->Point(ch,p2,AliHMPIDParam::kPc); //point & norm for entrance to PC plane |
97a4d538 | 100 | if(pTrk->Intersect(p1,n1)==kFALSE) return -1; //try to intersect track with the middle of radiator |
101 | if(pTrk->Intersect(p2,n2)==kFALSE) return -1; | |
496c71b0 | 102 | pParam->Mars2LorsVec(ch,n1,theta,phi); //track angles at RAD |
103 | pParam->Mars2Lors (ch,p1,xRa,yRa); //TRKxRAD position | |
104 | pParam->Mars2Lors (ch,p2,xPc,yPc); //TRKxPC position | |
105 | if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return ch; //return intersected chamber | |
106 | return -1; //no intersection with HMPID chambers | |
107 | }//IntTrkCha() | |
108 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
d3da6dc4 | 109 | Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree) |
110 | { | |
59d9d4b3 | 111 | // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event. |
d3da6dc4 | 112 | // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR() |
113 | // Returns: error code (currently ignored in AliReconstruction::RunTraking()) | |
ae5a42aa | 114 | for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i])); |
94b1fbfa | 115 | pCluTree->GetEntry(0); |
94b1fbfa | 116 | return 0; |
59d9d4b3 | 117 | }//LoadClusters() |
d3da6dc4 | 118 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
af885e0f | 119 | Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd) |
abb5f786 | 120 | { |
59d9d4b3 | 121 | // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event |
abb5f786 | 122 | // Agruments: pEsd - pointer to ESD |
123 | // Returns: error code | |
f455af6e | 124 | AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean |
49881df7 | 125 | AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1 |
0a7c7084 | 126 | if(!pNmeanEnt) AliError("No Nmean C6F14 "); |
127 | if(!pQthreEnt) AliError("No Qthre"); | |
94b1fbfa | 128 | |
afe12692 | 129 | return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject()); |
611e810d | 130 | }//PropagateBack() |
abb5f786 | 131 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
afe12692 | 132 | Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre) |
d3da6dc4 | 133 | { |
aed5907d | 134 | // Static method to reconstruct the Cherenkov angle for all valid tracks of a given event. |
59d9d4b3 | 135 | // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time) |
aed5907d | 136 | // Returns: error code, 0 if no errors |
137 | // | |
138 | // Algortihm: Loop over tracks | |
139 | // | |
140 | // 1. Find the closest MIP cluster using fast tracks extrapolation method | |
141 | // 2. Propagate track to the MIP cluster using the STEER method | |
142 | // 3. Update the track information with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation) | |
143 | // 4. Propagate back the constrained track to the radiator radius ( not exact yet) | |
144 | // 5. Propagation in the last 10 cm with the fast method | |
145 | // 6. Set ESDtrack information | |
146 | // 7. Calculate the Cherenkov angle | |
147 | ||
148 | const Double_t kMaxSnp=0.9; //maximal snp for prolongation | |
149 | const Double_t kdRadiator=10; // distance between radiator and the plane | |
496c71b0 | 150 | AliHMPIDRecon recon; //instance of reconstruction class, nothing important in ctor |
b8541917 | 151 | AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam |
39cd22e6 | 152 | Float_t xPc,yPc,xRa,yRa,theta,phi; |
aed5907d | 153 | |
43e51ec6 | 154 | Double_t cluLORS[2]={0}; |
496c71b0 | 155 | Int_t nMipClusTot=0; |
aed5907d | 156 | |
84770d7b | 157 | Double_t qthre = 0; Double_t nmean=0; Int_t hvsec=0; |
25692ca4 | 158 | // Int_t nClusCh[AliHMPIDParam::kMaxCh+1]; |
3e630437 | 159 | |
15f675cb | 160 | Bool_t tsRight = kTRUE; |
161 | ||
3e630437 | 162 | UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); // |
163 | UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); // | |
164 | UInt_t ts = pEsd->GetTimeStamp(); | |
496c71b0 | 165 | |
3e630437 | 166 | if(ts<tsmin || ts>tsmax) { |
167 | AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts)); | |
15f675cb | 168 | tsRight = kFALSE; |
3e630437 | 169 | } |
170 | ||
496c71b0 | 171 | for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event |
c38d443f | 172 | |
aed5907d | 173 | // Double_t bestChi2=99999;chi2=99999; //init. track matching params |
25692ca4 | 174 | // Double_t dmin=999999,bz=0,distCut=1,distParams[5]={1}; |
175 | Double_t dmin=999999,distCut=1,distParams[5]={1}; | |
43e51ec6 | 176 | |
177 | Bool_t isOkDcut=kFALSE; | |
178 | Bool_t isOkQcut=kFALSE; | |
179 | Bool_t isMatched=kFALSE; | |
180 | ||
496c71b0 | 181 | AliHMPIDCluster *bestHmpCluster=0x0; //the best matching cluster |
c38d443f | 182 | AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //get reconstructed track |
183 | ||
184 | if(!pTrk->IsOn(AliESDtrack::kTPCout)) continue; | |
185 | ||
186 | if(pTrk->IsOn(AliESDtrack::kTPCrefit)) continue; | |
79439569 | 187 | AliESDfriendTrack *ftrack= (AliESDfriendTrack *)pTrk->GetFriendTrack(); |
188 | if (!ftrack) continue; | |
189 | if (!ftrack->GetTPCOut()) continue; | |
496c71b0 | 190 | AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching |
aed5907d | 191 | AliHMPIDtrack *hmpTrkConstrained = 0; //create a hmpid track to be used for propagation and matching |
79439569 | 192 | hmpTrk->Set(ftrack->GetTPCOut()->GetX(), ftrack->GetTPCOut()->GetAlpha(),ftrack->GetTPCOut()->GetParameter(), ftrack->GetTPCOut()->GetCovariance()); |
193 | // | |
25692ca4 | 194 | //bz=AliTracker::GetBz(); |
aed5907d | 195 | |
196 | //initial flags for HMPID ESD infos | |
43e51ec6 | 197 | pTrk->SetHMPIDtrk(0,0,0,0); //no intersection found |
198 | pTrk->SetHMPIDmip(0,0,0,0); //store mip info in any case | |
199 | pTrk->SetHMPIDcluIdx(99,99999); //chamber not found, mip not yet considered | |
200 | pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed); //ring reconstruction not yet performed | |
496c71b0 | 201 | |
af291e40 | 202 | Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi); //find the intersected chamber for this track |
4ccd1903 | 203 | if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;} //no intersection at all, go after next track |
43e51ec6 | 204 | |
97a4d538 | 205 | pTrk->SetHMPIDtrk(xPc,yPc,theta,phi); //store initial infos |
43e51ec6 | 206 | pTrk->SetHMPIDcluIdx(ipCh,9999); //set chamber, index of cluster + cluster size |
496c71b0 | 207 | |
aed5907d | 208 | // track intersects the chamber ipCh: find the MIP |
496c71b0 | 209 | |
210 | TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters | |
211 | nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber | |
25692ca4 | 212 | // nClusCh[ipCh] = nMipClusTot; |
496c71b0 | 213 | |
4ccd1903 | 214 | if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;} |
43e51ec6 | 215 | |
aed5907d | 216 | Int_t index=-1; //index of the "best" matching cluster |
217 | // | |
218 | // 1. Find the closest MIP cluster using fast tracks extrapolation method | |
219 | // | |
220 | for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop | |
496c71b0 | 221 | |
aed5907d | 222 | AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster |
223 | // evaluate qThre | |
15f675cb | 224 | if(tsRight){ |
aed5907d | 225 | if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) { |
226 | qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts); // | |
227 | } else { // in the past just 1 qthre | |
228 | hvsec = pParam->InHVSector(pClu->Y()); // per chamber | |
229 | if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts); // | |
15f675cb | 230 | } |
231 | } else qthre = pParam->QCut(); | |
aed5907d | 232 | |
233 | // | |
496c71b0 | 234 | if(pClu->Q()<qthre) continue; //charge compartible with MIP clusters |
235 | isOkQcut = kTRUE; | |
aed5907d | 236 | // |
237 | cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y(); //get the LORS coordinates of the cluster | |
43e51ec6 | 238 | Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1])); |
239 | ||
240 | if(dist<dmin) { | |
241 | dmin = dist; | |
242 | index=iClu; | |
243 | bestHmpCluster=pClu; | |
244 | } | |
4ccd1903 | 245 | } // clusters loop |
b8541917 | 246 | |
aed5907d | 247 | // moved down |
248 | /*if(!isOkQcut) { | |
496c71b0 | 249 | pTrk->SetHMPIDsignal(pParam->kMipQdcCut); |
4ccd1903 | 250 | delete hmpTrk;hmpTrk=0x0; continue; |
aed5907d | 251 | }*/ |
252 | ||
253 | // | |
254 | // 2. Propagate track to the MIP cluster using the STEER method | |
79439569 | 255 | // |
95d2d29e | 256 | |
257 | if(!bestHmpCluster) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;} | |
258 | ||
79439569 | 259 | TVector3 vG = pParam->Lors2Mars(ipCh,bestHmpCluster->X(),bestHmpCluster->Y()); |
260 | Double_t gx = vG[0]; | |
261 | Double_t gy = vG[1]; | |
262 | Double_t gz = vG[2]; | |
263 | Double_t alpha=TMath::ATan2(gy,gx); | |
264 | Double_t radiusH=TMath::Sqrt(gy*gy+gx*gx); | |
aed5907d | 265 | if (!(hmpTrk->Rotate(alpha,kTRUE))) continue; |
266 | if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrk,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;} | |
79439569 | 267 | // |
aed5907d | 268 | // 3. Update the track with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation) |
269 | // | |
270 | AliExternalTrackParam trackC(*hmpTrk); | |
271 | Double_t posC[2]={0,gz}; | |
272 | Double_t covC[3]={0.1*0.1, 0, 0.1*0.1}; | |
273 | trackC.Update(posC,covC); | |
274 | // | |
275 | // 4. Propagate back the constrained track to the radiator radius ( not exact yet) | |
276 | // | |
277 | hmpTrkConstrained = new AliHMPIDtrack(*pTrk); | |
278 | hmpTrkConstrained->Set(trackC.GetX(), trackC.GetAlpha(),trackC.GetParameter(), trackC.GetCovariance()); | |
279 | if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrkConstrained,radiusH-kdRadiator,pTrk->GetMass(),1,kFALSE,kMaxSnp,1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;} | |
280 | // | |
281 | // 5. Propagation in the last 10 cm with the fast method | |
282 | // | |
283 | Float_t xPc0=0, yPc0=0; | |
284 | ||
7243a52d | 285 | IntTrkCha(ipCh, hmpTrk, xPc0,yPc0,xRa,yRa,theta,phi); |
286 | IntTrkCha(ipCh, hmpTrkConstrained, xPc,yPc,xRa,yRa,theta,phi); | |
aed5907d | 287 | // |
288 | // 6. Set ESDtrack information | |
289 | // | |
290 | Int_t cluSiz = bestHmpCluster->Size(); | |
291 | pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0); //store mip info in any case | |
292 | pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set chamber, index of cluster + cluster size | |
293 | pTrk->SetHMPIDtrk(xPc0,yPc0,theta,phi); | |
79439569 | 294 | // |
79439569 | 295 | // |
296 | // Dump debug info if specified | |
297 | // | |
aed5907d | 298 | if (AliHMPIDReconstructor::StreamLevel()>0) { |
79439569 | 299 | AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut())); |
300 | AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk); | |
bf3b953e | 301 | if(!trackTPC->Rotate(alpha)) continue; |
302 | if(!trackCurrent->Rotate(alpha)) continue; | |
79439569 | 303 | Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1); |
f23adc8f | 304 | Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1); |
305 | AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut())); | |
bf3b953e | 306 | if(!trackTPCNB->Rotate(alpha)) continue; |
f23adc8f | 307 | Bool_t statusTPCNB=kTRUE; |
308 | Double_t bfield[3]={0,0,0}; | |
309 | for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){ | |
310 | Double_t xyz[3]; | |
311 | trackTPCNB->GetXYZ(xyz); | |
312 | GetBxByBz(xyz,bfield); | |
313 | statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield); | |
314 | } | |
315 | statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield); | |
316 | // | |
aed5907d | 317 | Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp())); |
79439569 | 318 | Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz()); |
319 | // | |
aed5907d | 320 | AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC); |
321 | Double_t pos[2]={0,gz}; | |
322 | Double_t cov[3]={0.1*0.1, 0, 0.1*0.1}; | |
323 | Double_t chi2C = trackTPCConstrained->GetPredictedChi2(pos,cov); | |
324 | trackTPCConstrained->Update(pos,cov); | |
79439569 | 325 | (*fDebugStreamer)<<"track"<< |
aed5907d | 326 | "rH="<<radiusH<< // radius of cluster |
327 | "angle="<<tanAlpha<< // tan of the local inlination angle | |
328 | "dC="<<deltaC<< // delta of the curvature | |
f23adc8f | 329 | "trackTPC.="<<trackTPC<< // TPC outer param extrapolated to the HMPID |
330 | "trackTPCNB.="<<trackTPCNB<< // TPC track prpagated with material budget correction | |
aed5907d | 331 | "chi2C="<<chi2C<< |
f23adc8f | 332 | "trackTPCC.="<<trackTPCConstrained<< // TPC outer param extrapolated to the HMPID constrained |
aed5907d | 333 | "trackCurrent.="<<trackCurrent<< // current track extrapolated to the HMPID |
334 | "sTPC="<<statusTPC<< // status for the current TPC track | |
335 | "sCurrent="<<statusCurrent<< // status for the current global track | |
336 | "cl.="<<bestHmpCluster<< // HMPID cluster | |
79439569 | 337 | // |
aed5907d | 338 | "t.="<<pTrk<< // curent esd track |
339 | "ft.="<<ftrack<< // friend track | |
340 | "hmpTrk.="<<hmpTrk<< // hmpid tracks as used in the following code | |
f23adc8f | 341 | "hmpTrkC.="<<hmpTrkConstrained<< // constrained hmpid tracks as used in the following code |
aed5907d | 342 | "gx="<<gx<< // global cluster position X |
343 | "gy="<<gy<< // Y | |
344 | "gz="<<gz<< // Z | |
79439569 | 345 | "\n"; |
aed5907d | 346 | } |
347 | // | |
348 | // | |
349 | // | |
350 | if(!isOkQcut) { | |
351 | pTrk->SetHMPIDsignal(pParam->kMipQdcCut); | |
352 | delete hmpTrk;hmpTrk=0x0; | |
353 | delete hmpTrkConstrained;hmpTrkConstrained=0x0; | |
354 | continue; | |
355 | } | |
79439569 | 356 | |
b8541917 | 357 | if(AliHMPIDReconstructor::GetRecoParam()) //retrieve distance cut |
358 | { | |
b8541917 | 359 | if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE) //distance cut is fixed number |
360 | { | |
361 | distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist(); | |
362 | } | |
363 | else | |
364 | { | |
13fb7604 | 365 | for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar); //prevision: distance cut is function of momentum |
b8541917 | 366 | distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise |
367 | } | |
368 | } | |
369 | else {distCut=pParam->DistCut();} | |
aed5907d | 370 | |
371 | //dmin recalculated | |
372 | ||
373 | dmin = TMath::Sqrt((xPc0-bestHmpCluster->X())*(xPc0-bestHmpCluster->X())+(yPc0-bestHmpCluster->Y())*(yPc0-bestHmpCluster->Y())); | |
374 | ||
b8541917 | 375 | if(dmin < distCut) { |
496c71b0 | 376 | isOkDcut = kTRUE; |
4ccd1903 | 377 | } |
aed5907d | 378 | //isOkDcut = kTRUE; // switch OFF cut |
b8541917 | 379 | |
496c71b0 | 380 | if(!isOkDcut) { |
496c71b0 | 381 | pTrk->SetHMPIDsignal(pParam->kMipDistCut); //closest cluster with enough charge is still too far from intersection |
382 | } | |
383 | ||
384 | if(isOkQcut*isOkDcut) isMatched = kTRUE; // MIP-Track matched !! | |
385 | ||
aed5907d | 386 | if(!isMatched) {delete hmpTrk;hmpTrk=0x0;delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;} // If matched continue... |
43e51ec6 | 387 | |
aed5907d | 388 | Bool_t isOk = kTRUE; |
95d2d29e | 389 | if(!isOk) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;} |
aed5907d | 390 | pTrk->SetOuterHmpParam(hmpTrkConstrained,AliESDtrack::kHMPIDout); |
c38d443f | 391 | |
392 | FillResiduals(hmpTrk,bestHmpCluster,kFALSE); | |
b8541917 | 393 | |
aed5907d | 394 | Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved |
496c71b0 | 395 | |
496c71b0 | 396 | //evaluate nMean |
15f675cb | 397 | if(tsRight){ |
398 | if(pNmean->GetEntries()==21) { //for backward compatibility | |
399 | nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts); //C6F14 Nmean for this chamber | |
400 | } else { | |
15f675cb | 401 | if(iRad < 0) { |
0a7c7084 | 402 | nmean = -1; |
15f675cb | 403 | } else { |
404 | Double_t tLow = ((TF1*)pNmean->At(6*ipCh+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber | |
405 | Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber | |
406 | Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y | |
407 | nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp | |
408 | } | |
409 | if(nmean < 0){ //track didn' t pass through the radiator | |
496c71b0 | 410 | pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag |
411 | pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set index of cluster | |
4ccd1903 | 412 | delete hmpTrk;hmpTrk=0x0; |
aed5907d | 413 | delete hmpTrkConstrained;hmpTrkConstrained=0x0; |
496c71b0 | 414 | continue; |
15f675cb | 415 | } |
aed5907d | 416 | } |
417 | } else nmean = pParam->MeanIdxRad(); | |
418 | // | |
419 | // 7. Calculate the Cherenkov angle | |
496c71b0 | 420 | // |
aed5907d | 421 | recon.SetImpPC(xPc0,yPc0); //store track impact to PC |
97a4d538 | 422 | recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa); //search for Cerenkov angle of this track |
aed5907d | 423 | |
424 | Double_t thetaCkov = pTrk->GetHMPIDsignal(); | |
425 | ||
426 | if (AliHMPIDReconstructor::StreamLevel()>0) { | |
427 | AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut())); | |
428 | AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk); | |
bf3b953e | 429 | if(!trackTPC->Rotate(alpha)) continue; |
430 | if(!trackCurrent->Rotate(alpha)) continue; | |
aed5907d | 431 | Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1); |
432 | Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1); | |
433 | Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp())); | |
434 | Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz()); | |
435 | // | |
f23adc8f | 436 | AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut())); |
bf3b953e | 437 | if(!trackTPCNB->Rotate(alpha)) continue; |
f23adc8f | 438 | Bool_t statusTPCNB=kTRUE; |
439 | Double_t bfield[3]={0,0,0}; | |
440 | for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){ | |
441 | Double_t xyz[3]; | |
442 | trackTPCNB->GetXYZ(xyz); | |
443 | GetBxByBz(xyz,bfield); | |
444 | statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield); | |
445 | } | |
446 | statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield); | |
447 | ||
aed5907d | 448 | AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC); |
449 | Double_t pos[2]={0,gz}; | |
450 | Double_t cov[3]={0.1*0.1, 0, 0.1*0.1}; | |
451 | Double_t chi2C = trackTPCConstrained->GetPredictedChi2(pos,cov); | |
452 | trackTPCConstrained->Update(pos,cov); | |
453 | (*fDebugStreamer)<<"track2"<< | |
454 | "rH="<<radiusH<< // radius of cluster | |
455 | "angle="<<tanAlpha<< // tan of the local inlination angle | |
456 | "dC="<<deltaC<< // delta of the curvature | |
f23adc8f | 457 | "trackTPC.="<<trackTPC<< // TPC outer param extrapolated to the HMPID |
458 | "trackTPCNB.="<<trackTPCNB<< // TPC outer param extrapolated to the HMPID | |
aed5907d | 459 | "chi2C="<<chi2C<< |
f23adc8f | 460 | "trackTPCC.="<<trackTPCConstrained<< // TPC outer param extrapolated to the HMPID constrained |
aed5907d | 461 | "trackCurrent.="<<trackCurrent<< // current track extrapolated to the HMPID |
462 | "sTPC="<<statusTPC<< // status for the current TPC track | |
463 | "sCurrent="<<statusCurrent<< // status for the current global track | |
464 | "cl.="<<bestHmpCluster<< // HMPID cluster | |
465 | // | |
466 | "t.="<<pTrk<< // curent esd track | |
467 | "ft.="<<ftrack<< // friend track | |
468 | "hmpTrk.="<<hmpTrk<< // hmpid tracks as used in the following code | |
f23adc8f | 469 | "hmpTrkC.="<<hmpTrkConstrained<< // constrained hmpid tracks as used in the following code |
aed5907d | 470 | "gx="<<gx<< // global cluster position X |
471 | "gy="<<gy<< // Y | |
472 | "gz="<<gz<< // Z | |
473 | "thetaCkov="<<thetaCkov<< // Cherenkov angle | |
474 | "\n"; | |
475 | } | |
476 | ||
477 | if(pTrk->GetHMPIDsignal()<0) { | |
478 | delete hmpTrk;hmpTrk=0x0; | |
479 | delete hmpTrkConstrained; hmpTrkConstrained=0x0; | |
480 | continue;} | |
235a21c2 | 481 | |
e56b695f | 482 | AliHMPIDPid pID; |
483 | Double_t prob[5]; | |
b7a2d22d | 484 | pID.FindPid(pTrk,nmean,5,prob); |
e56b695f | 485 | pTrk->SetHMPIDpid(prob); |
aed5907d | 486 | delete hmpTrk; hmpTrk=0x0; |
487 | delete hmpTrkConstrained; hmpTrkConstrained=0x0; | |
496c71b0 | 488 | }//iTrk |
489 | ||
d3da6dc4 | 490 | return 0; // error code: 0=no error; |
611e810d | 491 | }//Recon() |
492 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
b7a2d22d | 493 | Int_t AliHMPIDTracker::ReconFromKin(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre) |
494 | { | |
495 | // Static method to reconstruct Theta Ckov for all valid tracks of a given event. | |
496 | // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time) | |
497 | // Returns: error code, 0 if no errors | |
498 | ||
499 | AliHMPIDRecon recon; //instance of reconstruction class, nothing important in ctor | |
500 | AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam | |
501 | Float_t xPc,yPc,xRa,yRa,theta,phi; | |
502 | Double_t cluLORS[2]={0}; | |
503 | Int_t nMipClusTot=0; | |
504 | Double_t qthre = 0; Double_t nmean=0; Int_t hvsec=0; | |
23874fc6 | 505 | //Int_t nClusCh[AliHMPIDParam::kMaxCh+1]; |
b7a2d22d | 506 | |
507 | Bool_t tsRight = kTRUE; | |
508 | ||
509 | UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); // | |
510 | UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); // | |
511 | UInt_t ts = pEsd->GetTimeStamp(); | |
512 | ||
513 | if(ts<tsmin || ts>tsmax) { | |
514 | AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts)); | |
515 | tsRight = kFALSE; | |
516 | } | |
517 | ||
518 | for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event | |
519 | ||
23874fc6 | 520 | Double_t dmin=999999,distCut=1,distParams[5]={1}; |
b7a2d22d | 521 | |
522 | Bool_t isOkDcut=kFALSE; | |
523 | Bool_t isOkQcut=kFALSE; | |
524 | Bool_t isMatched=kFALSE; | |
525 | ||
526 | AliHMPIDCluster *bestHmpCluster=0x0; //the best matching cluster | |
527 | AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //get reconstructed track | |
528 | ||
529 | AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching | |
23874fc6 | 530 | // bz=AliTracker::GetBz(); |
b7a2d22d | 531 | //initial flags for HMPID ESD infos |
532 | pTrk->SetHMPIDtrk(0,0,0,0); //no intersection found | |
533 | pTrk->SetHMPIDmip(0,0,0,0); //store mip info in any case | |
534 | pTrk->SetHMPIDcluIdx(99,99999); //chamber not found, mip not yet considered | |
535 | pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed); //ring reconstruction not yet performed | |
536 | ||
537 | Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi); //find the intersected chamber for this track | |
538 | if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;} //no intersection at all, go after next track | |
539 | ||
540 | pTrk->SetHMPIDtrk(xPc,yPc,theta,phi); //store initial infos | |
541 | pTrk->SetHMPIDcluIdx(ipCh,9999); //set chamber, index of cluster + cluster size | |
542 | ||
543 | // track intersects the chamber ipCh: find the MIP | |
544 | ||
545 | TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters | |
546 | nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber | |
23874fc6 | 547 | // nClusCh[ipCh] = nMipClusTot; |
b7a2d22d | 548 | |
549 | if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;} | |
550 | ||
551 | Int_t index=-1; //index of the "best" matching cluster | |
552 | ||
553 | for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop | |
554 | ||
555 | AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster | |
556 | // evaluate qThre | |
557 | if(tsRight){ | |
558 | if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) { | |
559 | qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts); // | |
560 | } else { // in the past just 1 qthre | |
561 | hvsec = pParam->InHVSector(pClu->Y()); // per chamber | |
562 | if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts); // | |
563 | } | |
564 | } else qthre = pParam->QCut(); | |
565 | ||
566 | // | |
567 | if(pClu->Q()<qthre) continue; //charge compartible with MIP clusters | |
568 | isOkQcut = kTRUE; | |
569 | // | |
570 | cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y(); //get the LORS coordinates of the cluster | |
571 | Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1])); | |
572 | ||
573 | if(dist<dmin) { | |
574 | dmin = dist; | |
575 | index=iClu; | |
576 | bestHmpCluster=pClu; | |
577 | } | |
578 | } // clusters loop | |
579 | ||
580 | if(!isOkQcut) { | |
581 | pTrk->SetHMPIDsignal(pParam->kMipQdcCut); | |
582 | delete hmpTrk;hmpTrk=0x0; continue; | |
583 | } | |
584 | ||
585 | Double_t radius = (pParam->Lors2Mars(ipCh,pParam->SizeAllX()/2,pParam->SizeAllY()/2)).Mag(); | |
586 | ||
587 | if(!AliTracker::PropagateTrackToBxByBz(hmpTrk,radius,pTrk->GetMass(),1,kFALSE)) {delete hmpTrk;hmpTrk=0x0;continue;} | |
588 | ||
589 | if(!hmpTrk->PropagateTo(bestHmpCluster)) {delete hmpTrk;hmpTrk=0x0;continue;} | |
590 | ||
591 | Int_t cluSiz = bestHmpCluster->Size(); | |
592 | pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0); //store mip info in any case | |
593 | pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set chamber, index of cluster + cluster size | |
594 | ||
595 | if(AliHMPIDReconstructor::GetRecoParam()) //retrieve distance cut | |
596 | { | |
597 | if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE) //distance cut is fixed number | |
598 | { | |
599 | distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist(); | |
600 | } | |
601 | else | |
602 | { | |
603 | for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar); //prevision: distance cut is function of momentum | |
604 | distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise | |
605 | } | |
606 | } | |
607 | else {distCut=pParam->DistCut();} | |
608 | ||
609 | if(dmin < distCut) { | |
610 | isOkDcut = kTRUE; | |
611 | } | |
612 | ||
613 | if(!isOkDcut) { | |
614 | pTrk->SetHMPIDsignal(pParam->kMipDistCut); //closest cluster with enough charge is still too far from intersection | |
615 | } | |
616 | ||
617 | if(isOkQcut*isOkDcut) isMatched = kTRUE; // MIP-Track matched !! | |
618 | ||
619 | if(!isMatched) {delete hmpTrk;hmpTrk=0x0;continue;} // If matched continue... | |
620 | ||
621 | Bool_t isOk = hmpTrk->Update(bestHmpCluster,0.1,0); | |
622 | if(!isOk) {delete hmpTrk;hmpTrk=0x0;continue;} | |
623 | pTrk->SetOuterHmpParam(hmpTrk,AliESDtrack::kHMPIDout); | |
624 | ||
625 | FillResiduals(hmpTrk,bestHmpCluster,kFALSE); | |
626 | ||
627 | Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved | |
628 | ||
629 | //evaluate nMean | |
630 | if(tsRight){ | |
631 | if(pNmean->GetEntries()==21) { //for backward compatibility | |
632 | nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts); //C6F14 Nmean for this chamber | |
633 | } else { | |
634 | if(iRad < 0) { | |
635 | nmean = -1; | |
636 | } else { | |
637 | Double_t tLow = ((TF1*)pNmean->At(6*ipCh+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber | |
638 | Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber | |
639 | Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y | |
640 | nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp | |
641 | } | |
642 | if(nmean < 0){ //track didn' t pass through the radiator | |
643 | pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag | |
644 | pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set index of cluster | |
645 | delete hmpTrk;hmpTrk=0x0; | |
646 | continue; | |
647 | } | |
648 | } | |
649 | } else nmean = pParam->MeanIdxRad(); | |
650 | ||
651 | // | |
652 | recon.SetImpPC(xPc,yPc); //store track impact to PC | |
653 | recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa); //search for Cerenkov angle of this track | |
654 | ||
655 | if(pTrk->GetHMPIDsignal()<0) {delete hmpTrk;hmpTrk=0x0;continue;} | |
656 | ||
657 | AliHMPIDPid pID; | |
658 | Double_t prob[5]; | |
659 | pID.FindPid(pTrk,nmean,5,prob); | |
660 | // Printf("Tracker RefIndex = %f, chmaber = %i, prob1 = %f, prob2 = %f, prob3 = %f, prob4 = %f, prob5 = %f", nmean,ipCh,prob[0],prob[1],prob[2],prob[3],prob[4]); | |
661 | pTrk->SetHMPIDpid(prob); | |
662 | delete hmpTrk;hmpTrk=0x0; | |
663 | }//iTrk | |
664 | ||
665 | return 0; // error code: 0=no error; | |
666 | }//ReconFromKin() | |
667 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
9785d5fb | 668 | Int_t AliHMPIDTracker::ReconHiddenTrk(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre) |
611e810d | 669 | { |
670 | // Static method to reconstruct Theta Ckov for all valid tracks of a given event. | |
afe12692 | 671 | // 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 | 672 | // Returns: error code, 0 if no errors |
9785d5fb | 673 | |
674 | AliHMPIDReconHTA reconHTA; //instance of reconstruction class, nothing important in ctor | |
675 | ||
676 | AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam | |
15f675cb | 677 | |
678 | Bool_t tsRight = kTRUE; | |
679 | ||
3e630437 | 680 | UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); // |
681 | UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); // | |
682 | UInt_t ts = pEsd->GetTimeStamp(); | |
683 | ||
684 | if(ts<tsmin || ts>tsmax) { | |
685 | AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts)); | |
15f675cb | 686 | tsRight = kFALSE; |
3e630437 | 687 | } |
688 | ||
9785d5fb | 689 | for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event |
690 | ||
691 | AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //here it is simulated or just empty track | |
692 | Int_t ipCh; | |
693 | ipCh = pTrk->GetHMPIDcluIdx();ipCh/=1000000; | |
694 | if(ipCh<0) continue; | |
695 | ||
696 | TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters | |
697 | Int_t nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber | |
698 | ||
699 | Double_t qMip=-1; | |
700 | Int_t chMip=-1; | |
ee2f3539 | 701 | Double_t xMip = 0; |
702 | Double_t yMip = 0; | |
703 | Int_t indMip = 0; | |
704 | Int_t cluMipSiz = 0; | |
9785d5fb | 705 | |
706 | for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop | |
707 | ||
708 | AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster | |
709 | Double_t qClus = pClu->Q(); | |
710 | if(qClus>qMip) { | |
711 | qMip = qClus; | |
712 | chMip = pClu->Ch(); | |
713 | xMip = pClu->X(); | |
714 | yMip = pClu->Y(); | |
715 | indMip = iClu; | |
716 | cluMipSiz = pClu->Size(); | |
717 | } | |
718 | }//clus loop | |
719 | ||
720 | if(chMip<0) return 1; | |
721 | ||
722 | Int_t hvsec; | |
ee2f3539 | 723 | Double_t qthre=0; |
15f675cb | 724 | |
9785d5fb | 725 | // evaluate qThre |
15f675cb | 726 | if(tsRight){ |
9785d5fb | 727 | if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) { // just for backward compatibility |
3e630437 | 728 | qthre=((TF1*)pQthre->At(chMip))->Eval(ts); // |
9785d5fb | 729 | } else { // in the past just 1 qthre |
730 | hvsec = pParam->InHVSector(yMip); // per chamber | |
3e630437 | 731 | if(hvsec>=0) qthre=((TF1*)pQthre->At(6*chMip+hvsec))->Eval(ts); // |
15f675cb | 732 | } |
733 | } else qthre = pParam->QCut(); | |
9785d5fb | 734 | // |
735 | if(qMip<qthre) { | |
3e630437 | 736 | pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case |
9785d5fb | 737 | pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz); |
738 | pTrk->SetHMPIDsignal(pParam->kMipQdcCut); | |
739 | return 1; //charge compatible with MIP clusters | |
740 | } | |
741 | ||
742 | pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case | |
743 | pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz); | |
744 | ||
745 | Double_t yRa = yMip; //just an approx... | |
746 | Double_t nmean; | |
ddb21a01 | 747 | |
748 | Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved | |
749 | ||
9785d5fb | 750 | //evaluate nMean |
15f675cb | 751 | if(tsRight){ |
9785d5fb | 752 | if(pNmean->GetEntries()==21) { //for backward compatibility |
3e630437 | 753 | nmean=((TF1*)pNmean->At(3*chMip))->Eval(ts); //C6F14 Nmean for this chamber |
9785d5fb | 754 | } else { |
9785d5fb | 755 | if(iRad < 0) { |
756 | nmean = -1; | |
757 | } else { | |
3e630437 | 758 | Double_t tLow = ((TF1*)pNmean->At(6*chMip+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber |
759 | Double_t tHigh = ((TF1*)pNmean->At(6*chMip+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber | |
9785d5fb | 760 | Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y |
761 | nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp | |
762 | } | |
763 | if(nmean < 0){ //track didn' t pass through the radiator | |
764 | pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag | |
765 | return 1; | |
15f675cb | 766 | } |
767 | } | |
768 | } else nmean = pParam->MeanIdxRad(); | |
9785d5fb | 769 | // |
ee2f3539 | 770 | if(!reconHTA.CkovHiddenTrk(pTrk,(TClonesArray *)pClus->At(ipCh),indMip,nmean)) { //search for track parameters and Cerenkov angle of this track |
771 | AliHMPIDPid pID; | |
772 | Double_t prob[5]; | |
b7a2d22d | 773 | pID.FindPid(pTrk,nmean,5,prob); |
ee2f3539 | 774 | pTrk->SetHMPIDpid(prob); |
775 | } | |
e56b695f | 776 | // 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 | 777 | }//iTrk |
778 | ||
779 | return 0; // error code: 0=no error; | |
780 | ||
611e810d | 781 | }//Recon() |
d3da6dc4 | 782 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
c1af14f7 | 783 | void AliHMPIDTracker::FillClusterArray(TObjArray* array) const { |
784 | ||
785 | // Publishes all pointers to clusters known to the tracker into the | |
786 | // passed object array. | |
787 | // The ownership is not transfered - the caller is not expected to delete | |
788 | // the clusters | |
789 | ||
790 | for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){ | |
791 | TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh]; | |
792 | for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){ | |
793 | AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu); | |
794 | array->AddLast(pClu); | |
795 | }//cluster loop in iCh | |
796 | pCluArr->Delete(); | |
797 | }//Ch loop | |
798 | ||
799 | return; | |
800 | } | |
801 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |