ATLAS Offline Software
Loading...
Searching...
No Matches
ALFA_HalfReco.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
8 AthMessaging("ALFA_HalfReco")
9{
10 ATH_MSG_DEBUG("begin ALFA_HalfReco::ALFA_HalfReco");
11
12 m_fOvU = 0.0;
13 m_fOvV = 0.0;
14 m_fOverlapCut = 0.0;
15 m_fRecXPos = -9999.0;
16 m_fRecYPos = -9999.0;
17 m_iHalf = 0;
19
20// memset(&m_iNumHitsLayer, 0, sizeof(m_iNumHitsLayer));
21
22 memset(&m_fhits, 0, sizeof(m_fhits));
23 std::fill_n(m_fhits, sizeof(m_fhits)/sizeof(Int_t), -9999);
24
25 m_iNumU = 0;
26 m_iNumV = 0;
27 m_iRPot = 0;
28 m_iUVCut = 0;
29
30 ATH_MSG_DEBUG("end ALFA_HalfReco::ALFA_HalfReco");
31}
32
36
37StatusCode ALFA_HalfReco::Initialize(Float_t faMD[RPOTSCNT][ALFALAYERSCNT*ALFAPLATESCNT][ALFAFIBERSCNT], Float_t fbMD[RPOTSCNT][ALFALAYERSCNT*ALFAPLATESCNT][ALFAFIBERSCNT], Int_t iHalf, Int_t fMultiplicityCut, Int_t iUVCut, Float_t fOverlapCut)
38{
39 ATH_MSG_DEBUG("begin ALFA_HalfReco::Initialize()");
40
41 m_iMultiplicityCut = fMultiplicityCut;
42 m_iUVCut = iUVCut;
43 m_fOverlapCut = fOverlapCut;
44 m_iHalf = iHalf;
45
46 for (Int_t iRPot = 0; iRPot < RPOTSCNT; iRPot++)
47 {
48 for (Int_t iLayer = 0; iLayer < ALFALAYERSCNT*ALFAPLATESCNT; iLayer++)
49 {
50 for (Int_t iFiber = 0; iFiber < ALFAFIBERSCNT; iFiber++)
51 {
52 m_faMD[iRPot][iLayer][iFiber] = faMD[iRPot][iLayer][iFiber];
53 m_fbMD[iRPot][iLayer][iFiber] = fbMD[iRPot][iLayer][iFiber];
54 }
55 }
56 }
57
58 ATH_MSG_DEBUG("end ALFA_HalfReco::Initialize()");
59
60 return StatusCode::SUCCESS;
61}
62
63StatusCode ALFA_HalfReco::Execute(Int_t iRPot, const std::list<MDHIT> &ListMDHits)
64{
65 ATH_MSG_DEBUG("ALFA_HalfReco::Execute()");
66
67 FIBERS structFibers;
68 Int_t iNumUFiberHits = 0; //u_hits, v_hits (a condition for the UVCut 3)
69 Int_t iNumVFiberHits = 0;
70
71 m_iRPot = iRPot;
72 m_MapLayers.clear();
73
74 std::list<MDHIT>::const_iterator iter;
75 for (iter=ListMDHits.begin(); iter!=ListMDHits.end(); ++iter)
76 {
77 if (iRPot == (*iter).iRPot)
78 {
79 if(m_MapLayers.find((*iter).iPlate)==m_MapLayers.end())
80 {
81 structFibers.ListFibers.clear();
82 m_MapLayers.insert(std::pair<int, FIBERS>((*iter).iPlate, structFibers));
83 m_MapLayers[(*iter).iPlate].ListFibers.push_back((*iter).iFiber);
84 }
85 else
86 {
87 m_MapLayers[(*iter).iPlate].ListFibers.push_back((*iter).iFiber);
88 }
89 }
90 }
91
92 //Checking that the multiplicity cut conditions are satisfied
93 //At least more than UV_cut layers have a multiplicity lower than multi_cut
94// for (Int_t i=half_ch*10;i<half_ch*10+10;i++)
96 {
97// if (m_iNumHitsLayer[i] <= m_iMultiplicityCut && m_iNumHitsLayer[i] > 0)
98 if ((Int_t)m_MapLayers[i].ListFibers.size()<=m_iMultiplicityCut && (Int_t)m_MapLayers[i].ListFibers.size()>0)
99 {
100 if (fmod(double(i),double(2)) == 1) iNumUFiberHits++;
101 else iNumVFiberHits++;
102 }
103 }
104
105
106 // RECONSTRUCTION
107 m_fRecXPos = -9999.0;
108 m_fRecYPos = -9999.0;
109 m_fOvU = -9999.0;
110 m_fOvV = -9999.0;
111
112 memset(&m_fhits, 0, sizeof(m_fhits));
113 std::fill_n(m_fhits, sizeof(m_fhits)/sizeof(Int_t), -9999);
114
115 if (iNumUFiberHits>=m_iUVCut && iNumVFiberHits>=m_iUVCut) OverLap();
116
117 ATH_MSG_DEBUG("end ALFA_HalfReco::Execute()");
118 return StatusCode::SUCCESS;
119}
120
121StatusCode ALFA_HalfReco::Finalize(Float_t &fRecXPos, Float_t &fRecYPos)
122{
123 ATH_MSG_DEBUG("begin ALFA_HalfReco::Finalize()");
124
125 fRecXPos = m_fRecXPos;
126 fRecYPos = m_fRecYPos;
127
128 ATH_MSG_DEBUG("end ALFA_HalfReco::Execute()");
129 return StatusCode::SUCCESS;
130}
131
132void ALFA_HalfReco::HistFill(Float_t &b_p, Float_t &b_n, Float_t &Ov_p, Float_t &Ov_n, Int_t &NumU, Int_t &NumV, Int_t iFlag)
133{
134 ATH_MSG_DEBUG("ALFA_HalfReco::HistFill()");
135
136 Int_t iHit;
137 Float_t fFibCen;
138 Float_t b_ref[2] = {-127.0, -127.0};
139 Float_t fXInter, fYInter;
140 Float_t p_min, p_max;
141 Float_t n_min, n_max;
142 Int_t iMinP_tmp, iMaxP_tmp;
143 Int_t iMinN_tmp, iMaxN_tmp;
144 Int_t iFullWidth;
145 const Int_t NBINTOT = 72000;
146 std::vector<Int_t> Over_p(NBINTOT);
147 std::vector<Int_t> Over_n(NBINTOT);
148
149 if (iFlag==0) //First call, Event not reconstructed for U and V
150 {
151 b_ref[0]=-127.0;
152 b_ref[1]=-127.0;
153 }
154
155 if (iFlag==1) //Event reconstructed for U but not for V so let's try to refine
156 {
157 b_ref[0]=b_p;
158 b_ref[1]=-127.0;
159
160 }
161
162 if (iFlag==2) //Event reconstructed for V but not for U so let's try to refine
163 {
164 b_ref[0]=-127.0;
165 b_ref[1]=b_n;
166 }
167
168 if (iFlag==3) //Event reconstructed for U and V, ultime refinement
169 {
170 b_ref[0]=b_p;
171 b_ref[1]=b_n;
172 }
173
174 for (Int_t i=0;i<72000;i++)
175 {
176 Over_p[i]=0;
177 Over_n[i]=0;
178 }
179
180 std::list<int>::iterator intIter;
181//for(Int_t i=half_ch*10;i<half_ch*10+10;i++)
182 for (Int_t iLayer = ALFAPLATESCNT*m_iHalf; iLayer < m_iHalf*ALFAPLATESCNT+ALFAPLATESCNT; iLayer++)
183 {
184 for (intIter=m_MapLayers[iLayer].ListFibers.begin(); intIter!=m_MapLayers[iLayer].ListFibers.end(); ++intIter)
185 {
186 iHit = *intIter;
187
188 if (m_faMD[m_iRPot][iLayer][iHit] >= 0)
189 {
190 fXInter = (b_ref[1]-m_fbMD[m_iRPot][iLayer][iHit])/(1+m_faMD[m_iRPot][iLayer][iHit]);
191 fYInter = (m_faMD[m_iRPot][iLayer][iHit]*b_ref[1]+m_fbMD[m_iRPot][iLayer][iHit])/(1+m_faMD[m_iRPot][iLayer][iHit])-b_ref[1];
192 fFibCen = (fYInter-fXInter)/sqrt(2.0);
193
194 for (Int_t bin=0;bin<480;bin++)
195 {
196 Over_p[(int)((fFibCen-0.24)*1000)+bin+36000]++;
197 }
198 }
199 else
200 {
201 fXInter = (m_fbMD[m_iRPot][iLayer][iHit]-b_ref[0])/(1-m_faMD[m_iRPot][iLayer][iHit]);
202 fYInter = (m_fbMD[m_iRPot][iLayer][iHit]-m_faMD[m_iRPot][iLayer][iHit]*b_ref[0])/(1-m_faMD[m_iRPot][iLayer][iHit])-b_ref[0];
203 fFibCen = (fYInter+fXInter)/sqrt(2.0);
204
205 for (Int_t bin=0;bin<480;bin++)
206 {
207 Over_n[(int)((fFibCen-0.24)*1000)+bin+36000]++;
208 }
209 }
210 }
211 }
212
213 //Determine maximum
214 NumU=0;
215 NumV=0;
216
217 std::vector<int> *max_p = new std::vector<int>;
218 std::vector<int> *max_n = new std::vector<int>;
219
220 for (Int_t i=0;i<72000;i++)
221 {
222 if (Over_p[i]>NumU) NumU = Over_p[i];
223 if (Over_n[i]>NumV) NumV = Over_n[i];
224 }
225
226 if (NumU >= m_iUVCut) // If more than m_iUVCut fibers in the maximum we do the reconstruction
227 {
228 for (Int_t i=0;i<72000;i++)
229 {
230 if (Over_p[i]==NumU)
231 {
232 max_p->push_back(i);
233 }
234 }
235
236 p_min = -36.0 + double(max_p->front())*1e-3;
237 p_max = -36.0 + double(max_p->back())*1e-3;
238 iMinP_tmp = max_p->front();
239 iMaxP_tmp = max_p->back();
240
241// Making sure that the plateau belongs to the same track
242 iFullWidth = max_p->size();
243
244 for (Int_t i=0; i<iFullWidth; i++)
245 {
246 if (max_p->at(iFullWidth-i-1)-iMinP_tmp < 500)
247 {
248 iMaxP_tmp = max_p->at(iFullWidth-i-1);
249 p_max = -36.0 + double(iMaxP_tmp)*1e-3;
250 break;
251 }
252 }
253
254 if ((p_max-p_min)<m_fOverlapCut) //Check for multiple tracks
255 {
256 b_p = b_ref[1]+(p_min+p_max)/sqrt(2.);
257 Ov_p = p_max-p_min;
258 }
259 else
260 {
261 b_p = -9999.0;
262 Ov_p = -9999.0;
263 }
264 }
265 else
266 {
267 b_p = -9999.0;
268 Ov_p = -9999.0;
269 }
270
271 if (NumV >= m_iUVCut) // If more than m_iUVCut fibers in the maximum we do the reconstruction
272 {
273 for (Int_t i=0;i<72000;i++)
274 {
275 if (Over_n[i]==NumV)
276 {
277 max_n->push_back(i);
278 }
279 }
280
281 n_min = -36.0 + double(max_n->front())*1e-3;
282 n_max = -36.0 + double(max_n->back())*1e-3;
283 iMinN_tmp = max_n->front();
284 iMaxN_tmp = max_n->back();
285
286// Making sure that the plateau belongs to the same track
287 iFullWidth = max_n->size();
288
289 for (Int_t i=0; i<iFullWidth; i++)
290 {
291 if (max_n->at(iFullWidth-i-1)-iMinN_tmp < 500)
292 {
293 iMaxN_tmp = max_n->at(iFullWidth-i-1);
294 n_max = -36.0 + double(iMaxN_tmp)*1e-3;
295 break;
296 }
297 }
298
299 if ((n_max-n_min)<m_fOverlapCut) //Check for multiple tracks
300 {
301 b_n = b_ref[0]+(n_min+n_max)/sqrt(2.);
302 Ov_n = n_max-n_min;
303 }
304 else
305 {
306 b_n = -9999.0;
307 Ov_n = -9999.0;
308 }
309 }
310 else
311 {
312 b_n = -9999.0;
313 Ov_n = -9999.0;
314 }
315
316 delete max_p;
317 delete max_n;
318
319 ATH_MSG_DEBUG("endl ALFA_HalfReco::HistFill()");
320}
321
323{
324 ATH_MSG_DEBUG("ALFA_HalfReco::OverLap()");
325
326 Float_t b_mean_n=-9999.0;
327 Float_t b_mean_p=-9999.0;
328 Float_t overReg_p;
329 Float_t overReg_n;
330 Int_t NumU, NumV;
331
332 //First call, projection on middle fiber for U and V
333 HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 0);
334
335 //If first projection worked, we do the ultimate one
336 if (b_mean_p!=-9999.0 && b_mean_n!=-9999.0)
337 {
338 HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 3);
339 if (b_mean_p!=-9999.0 && b_mean_n!=-9999.0) HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 3);
340 }
341 else
342 {
343 //If first projection worked for p side but not for n, we try to refine n projection using the information from p side
344 if (b_mean_p!=-9999.0 && b_mean_n==-9999.0)
345 {
346 HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 1);
347
348 if (b_mean_p!=-9999.0 && b_mean_n!=-9999.0)
349 {
350 //if it as worked for both we do the ultimate one
351 HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 3);
352 if (b_mean_p!=-9999.0 && b_mean_n!=-9999.0) HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 3);
353 }
354 }
355
356 //If first projection worked for p side but not for n, we try to refine n projection using the information from p side
357 if (b_mean_p==-9999.0 && b_mean_n!=-9999.0)
358 {
359 HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 2);
360
361 if (b_mean_p!=-9999.0 && b_mean_n!=-9999.0)
362 {
363 //if it as worked for both we do the ultimate one
364 HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 3);
365 if (b_mean_p!=-9999.0 && b_mean_n!=-9999.0) HistFill(b_mean_p, b_mean_n, overReg_p, overReg_n, NumU, NumV, 3);
366 }
367 }
368 }
369
370 Float_t x_tmp, y_tmp;
371 Float_t min_x,min_y;
372 Float_t dist_x,dist_y;
373// Float_t dist_min_x,dist_min_y;
374 Int_t gFib_x, gFib_y;
375 Int_t iHit;
376
377 Float_t OL_U, OL_V;
378
379 std::list<int>::iterator intIter;
380 if (b_mean_p!=-9999.0 && b_mean_n!=-9999.0)
381 {
382 m_fRecXPos = (b_mean_n-b_mean_p)/2.0;
383 m_fRecYPos = (b_mean_n+b_mean_p)/2.0;
384 OL_U = overReg_p;
385 OL_V = overReg_n;
386
387 //for (Int_t i=half_ch*10;i<half_ch*10+10;i++)
388 for (Int_t iLayer = ALFAPLATESCNT*m_iHalf; iLayer<m_iHalf*ALFAPLATESCNT+ALFAPLATESCNT; iLayer++)
389 {
390 min_x = 0.24/cos(TMath::Pi()/4.);
391 min_y = 0.24/cos(TMath::Pi()/4.);
392
393 gFib_x = 9999;
394 gFib_y = 9999;
395
396 for (intIter=m_MapLayers[iLayer].ListFibers.begin(); intIter!=m_MapLayers[iLayer].ListFibers.end(); ++intIter)
397 {
398 iHit = *intIter;
399 x_tmp = (m_fRecYPos-m_fbMD[m_iRPot][iLayer][iHit])/m_faMD[m_iRPot][iLayer][iHit];
400 y_tmp = m_faMD[m_iRPot][iLayer][iHit]*m_fRecXPos+m_fbMD[m_iRPot][iLayer][iHit];
401 dist_x = TMath::Abs(m_fRecXPos-x_tmp);
402 dist_y = TMath::Abs(m_fRecYPos-y_tmp);
403
404 if (dist_x <= min_x)
405 {
406// dist_min_x = m_fRecXPos-x_tmp;
407 min_x = dist_x;
408 gFib_x = iHit;
409 }
410
411 if (dist_y <= min_y)
412 {
413// dist_min_y = m_fRecYPos-y_tmp;
414 min_y = dist_y;
415 gFib_y = iHit;
416 }
417 }
418
419 if (gFib_x==gFib_y) m_fhits[iLayer]=gFib_x;
420 }
421 }
422 else
423 {
424 m_fRecXPos = -9999.0;
425 m_fRecYPos = -9999.0;
426 NumU = 0;
427 NumV = 0;
428 OL_U = -9999.0;
429 OL_V = -9999.0;
430 }
431
432 SetData(NumU, NumV, OL_U, OL_V);
433
434 ATH_MSG_DEBUG("end ALFA_HalfReco::OverLap()");
435}
436
437void ALFA_HalfReco::SetData(Int_t iNumU, Int_t iNumV, Float_t fOvU, Float_t fOvV)
438{
439 ATH_MSG_DEBUG("begin ALFA_HalfReco::SetData()");
440
441 m_iNumU = iNumU;
442 m_iNumV = iNumV;
443 m_fOvU = fOvU;
444 m_fOvV = fOvV;
445
446 ATH_MSG_DEBUG("end ALFA_HalfReco::SetData()");
447}
448
449void ALFA_HalfReco::GetData(Int_t &iNumU, Int_t &iNumV, Float_t &fOvU, Float_t &fOvV, Int_t (&iFibSel)[ALFALAYERSCNT*ALFAPLATESCNT])
450{
451 ATH_MSG_DEBUG("begin ALFA_HalfReco::GetData()");
452
453 iNumU = m_iNumU;
454 iNumV = m_iNumV;
455 fOvU = m_fOvU;
456 fOvV = m_fOvV;
457
458 for (Int_t iLayer=0; iLayer<ALFALAYERSCNT*ALFAPLATESCNT; iLayer++)
459 {
460 iFibSel[iLayer] = m_fhits[iLayer];
461 }
462
463 ATH_MSG_DEBUG("end ALFA_HalfReco::GetData()");
464}
#define RPOTSCNT
struct _FIBERS FIBERS
#define ALFALAYERSCNT
#define ALFAPLATESCNT
#define ALFAFIBERSCNT
#define ATH_MSG_DEBUG(x)
void SetData(Int_t iNumU, Int_t iNumV, Float_t fOvU, Float_t fOvV)
StatusCode Finalize(Float_t &fRecXPos, Float_t &fRecYPos)
Int_t m_iMultiplicityCut
Float_t m_fbMD[RPOTSCNT][ALFALAYERSCNT *ALFAPLATESCNT][ALFAFIBERSCNT]
StatusCode Execute(Int_t iRPot, const std::list< MDHIT > &ListMDHits)
StatusCode Initialize(Float_t faMD[RPOTSCNT][ALFALAYERSCNT *ALFAPLATESCNT][ALFAFIBERSCNT], Float_t fbMD[RPOTSCNT][ALFALAYERSCNT *ALFAPLATESCNT][ALFAFIBERSCNT], Int_t iHalf, Int_t fMultiplicityCut, Int_t iUVCut, Float_t fOverlapCut)
void GetData(Int_t &iNumU, Int_t &iNumV, Float_t &fOvU, Float_t &fOvV, Int_t(&iFibSel)[ALFALAYERSCNT *ALFAPLATESCNT])
Float_t m_fRecYPos
Float_t m_fRecXPos
Float_t m_fOverlapCut
void HistFill(Float_t &b_p, Float_t &b_n, Float_t &Ov_p, Float_t &Ov_n, Int_t &NumU, Int_t &NumV, Int_t iFlag)
Float_t m_faMD[RPOTSCNT][ALFALAYERSCNT *ALFAPLATESCNT][ALFAFIBERSCNT]
Int_t m_fhits[ALFALAYERSCNT *ALFAPLATESCNT]
std::map< int, FIBERS > m_MapLayers
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
std::list< int > ListFibers