EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CbmRichRingFinderHoughImpl.cxx
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file CbmRichRingFinderHoughImpl.cxx
1 
9 #include "CbmRichRingLight.h"
10 
11 #include "../../littrack/std/utils/CbmLitMemoryManagment.h"
12 #include "CbmRichRingFitterCOP.h"
13 #include "CbmRichRingSelectAnn.h"
14 
15 #include "TSystem.h"
16 
17 #include <iostream>
18 #include <cmath>
19 #include <set>
20 #include <algorithm>
21 #include <iostream>
22 
23 using std::cout;
24 using std::endl;
25 using std::vector;
26 using std::set;
27 using std::sort;
28 
30  fNofParts(0),
31 
32  fMaxDistance(0.f),
33  fMinDistance(0.f),
34  fMinDistanceSq(0.f),
35  fMaxDistanceSq(0.f),
36 
37  fMinRadius(0.f),
38  fMaxRadius(0.f),
39 
40  fDx(0.f),
41  fDy(0.f),
42  fDr(0.f),
43  fNofBinsX(0),
44  fNofBinsY(0),
45  fNofBinsXY(0),
46 
47  fHTCut(0),
48 
49  fNofBinsR(0),
50  fHTCutR(0),
51 
52  fMinNofHitsInArea(0),
53 
54  fRmsCoeffEl(0.f),
55  fMaxCutEl(0.f),
56  fRmsCoeffCOP(0.f),
57  fMaxCutCOP(0.f),
58 
59  fAnnCut(0.f),
60  fUsedHitsAllCut(0.f),
61 
62  fCurMinX(0.f),
63  fCurMinY(0.f),
64 
65  fData(),
66  fHist(),
67  fHistR(),
68  fHitInd(),
69  fFoundRings(),
70  fFitCOP(NULL),
71  fANNSelect(NULL)
72 {
73 
74 }
75 
77 {
78  if (NULL != fFitCOP) delete fFitCOP;
79  if (NULL != fANNSelect) delete fANNSelect;
80 }
81 
83 {
84  SetParameters();
85 
86  fHist.resize(fNofBinsXY);
87  fHistR.resize(fNofBinsR);
88  fHitInd.resize(fNofParts);
89 
92  fANNSelect->Init();
93 }
94 
96 {
97  if (fData.size() > MAX_NOF_HITS) {
98  cout << endl << endl << "-E- CbmRichRingFinderHoughImpl::DoFind"
99  << ". Number of hits is more than " << MAX_NOF_HITS << endl << endl;
100  return;
101  }
102 
103  for_each(fFoundRings.begin(), fFoundRings.end(), DeleteObject());
104  fFoundRings.clear();
105  fFoundRings.reserve(100);
106 
107  std::sort(fData.begin(), fData.end(), CbmRichHoughHitCmpUp());
109  RingSelection();
110  fData.clear();
111 }
112 
114 {
115  fMaxDistance = 11.5;
116  fMinDistance = 3.;
119 
120  fMinRadius = 3.3;
121  fMaxRadius = 5.7;
122 
123  fHTCut = 7;
124  fHTCutR = 5;
125  fMinNofHitsInArea = 4;
126 
127  fNofBinsX = 25;
128  fNofBinsY = 25;
129  fNofBinsR = 32;
130 
131  fAnnCut = 0.6;
132  fUsedHitsAllCut = 0.4;
133 
134  fRmsCoeffEl = 2.5;
135  fMaxCutEl = 1.0;
136  fRmsCoeffCOP = 3.;
137  fMaxCutCOP = 1.0;
138 
139  fNofParts = 2;
140  fDx = 2.f * fMaxDistance / (float)fNofBinsX;
141  fDy = 2.f * fMaxDistance / (float)fNofBinsY;
142  fDr = fMaxRadius / (float)fNofBinsR;
144 }
145 
147 {
148  int indmin, indmax;
149  unsigned int size = fData.size();
150  for (unsigned int iHit = 0; iHit < size; iHit++){
151  if (fData[iHit].fIsUsed == true) continue;
152 
153  fCurMinX = fData[iHit].fHit.fX - fMaxDistance;
154  fCurMinY = fData[iHit].fHit.fY - fMaxDistance;
155 
156  DefineLocalAreaAndHits(fData[iHit].fHit.fX, fData[iHit].fHit.fY , &indmin, &indmax);
157  HoughTransform(indmin, indmax);
158  FindPeak(indmin, indmax);
159  }
160 }
161 
163  float x0,
164  float y0,
165  int *indmin,
166  int *indmax)
167 {
168  CbmRichHoughHit mpnt;
169  vector<CbmRichHoughHit>::iterator itmin, itmax;
170 
171  //find all hits which are in the corridor
172  mpnt.fHit.fX = x0 - fMaxDistance;
173  itmin = std::lower_bound(fData.begin(), fData.end(), mpnt, CbmRichHoughHitCmpUp());
174 
175  mpnt.fHit.fX = x0 + fMaxDistance;
176  itmax = std::lower_bound(fData.begin(), fData.end(), mpnt, CbmRichHoughHitCmpUp()) - 1;
177 
178  *indmin = itmin - fData.begin();
179  *indmax = itmax - fData.begin();
180 
181  int arSize = *indmax - *indmin + 1;
182  if (arSize <= fMinNofHitsInArea) return;
183 
184  for (unsigned short i = 0; i < fNofParts; i++){
185  fHitInd[i].clear();
186  fHitInd[i].reserve( (*indmax-*indmin) / fNofParts);
187  }
188 
189  unsigned short indmin1 = (unsigned short)(*indmin);
190  unsigned short indmax1 = (unsigned short)(*indmax);
191 
192  for (unsigned short i = indmin1; i <= indmax1; i++) {
193  if (fData[i].fIsUsed == true) continue;
194  float ry = y0 - fData[i].fHit.fY;
195  if (fabs(ry) > fMaxDistance) continue;
196  float rx = x0 - fData[i].fHit.fX;
197  float d = rx * rx + ry * ry;
198  if (d > fMaxDistanceSq) continue;
199  fHitInd[i % fNofParts].push_back(i);
200  }
201 
202  for (unsigned short j = 0; j < fNofBinsXY; j++){
203  fHist[j] = 0;
204  }
205 
206  for (unsigned short k = 0; k < fNofBinsR; k++) {
207  fHistR[k] = 0;
208  }
209 }
210 
212  unsigned short int indmin,
213  unsigned short int indmax)
214 {
215  for (int iPart = 0; iPart < fNofParts; iPart++){
216  HoughTransformGroup(indmin, indmax, iPart);
217  }//iPart
218 }
219 
221  unsigned short int indmin,
222  unsigned short int indmax,
223  int iPart)
224 {
225  unsigned short nofHits = fHitInd[iPart].size();
226  float xcs, ycs; // xcs = xc - fCurMinX
227  float dx = 1.0f/fDx, dy = 1.0f/fDy, dr = 1.0f/fDr;
228 
229  vector<CbmRichHoughHit> data;
230  data.resize(nofHits);
231  for (int i = 0; i < nofHits; i++){
232  data[i] = fData[ fHitInd[iPart][i] ];
233  }
234 
235  typedef vector<CbmRichHoughHit>::iterator iH;
236 
237  for (iH iHit1 = data.begin(); iHit1 != data.end(); iHit1++) {
238  float iH1X = iHit1->fHit.fX;
239  float iH1Y = iHit1->fHit.fY;
240 
241  for (iH iHit2 = iHit1 + 1; iHit2 != data.end(); iHit2++) {
242  float iH2X = iHit2->fHit.fX;
243  float iH2Y = iHit2->fHit.fY;
244 
245  float rx0 = iH1X - iH2X;//rx12
246  float ry0 = iH1Y- iH2Y;//ry12
247  float r12 = rx0 * rx0 + ry0 * ry0;
248  if (r12 < fMinDistanceSq || r12 > fMaxDistanceSq) continue;
249 
250  float t10 = iHit1->fX2plusY2 - iHit2->fX2plusY2;
251  for (iH iHit3 = iHit2 + 1; iHit3 != data.end(); iHit3++) {
252  float iH3X = iHit3->fHit.fX;
253  float iH3Y = iHit3->fHit.fY;
254 
255  float rx1 = iH1X - iH3X;//rx13
256  float ry1 = iH1Y - iH3Y;//ry13
257  float r13 = rx1 * rx1 + ry1 * ry1;
258  if (r13 < fMinDistanceSq || r13 > fMaxDistanceSq)continue;
259 
260  float rx2 = iH2X - iH3X;//rx23
261  float ry2 = iH2Y - iH3Y;//ry23
262  float r23 = rx2 * rx2 + ry2 * ry2;
263  if (r23 < fMinDistanceSq || r23 > fMaxDistanceSq)continue;
264 
265  float det = rx2*ry0 - rx0*ry2;
266  if (det == 0.0f) continue;
267  float t19 = 0.5f / det;
268  float t5 = iHit2->fX2plusY2 - iHit3->fX2plusY2;
269 
270  float xc = (t5 * ry0 - t10 * ry2) * t19;
271  xcs = xc - fCurMinX;
272  int intX = int( xcs *dx);
273  if (intX < 0 || intX >= fNofBinsX ) continue;
274 
275  float yc = (t10 * rx2 - t5 * rx0) * t19;
276  ycs = yc - fCurMinY;
277  int intY = int( ycs *dy);
278  if (intY < 0 || intY >= fNofBinsY ) continue;
279 
280  //radius calculation
281  float t6 = iH1X - xc;
282  float t7 = iH1Y - yc;
283  //if (t6 > fMaxRadius || t7 > fMaxRadius) continue;
284  float r = sqrt(t6 * t6 + t7 * t7);
285  //if (r < fMinRadius) continue;
286  int intR = int(r *dr);
287  if (intR < 0 || intR >= fNofBinsR) continue;
288 
289  fHist[intX * fNofBinsX + intY]++;
290  fHistR[intR]++;
291  }// iHit1
292  }// iHit2
293  }// iHit3
294  //hitIndPart.clear();
295 }
296 
298  int indmin,
299  int indmax)
300 {
301  // Find MAX bin R and compare it with CUT
302  int maxBinR = -1, maxR = -1;
303  unsigned int size = fHistR.size();
304  for (unsigned int k = 0; k < size; k++){
305  if (fHistR[k] > maxBinR){
306  maxBinR = fHistR[k];
307  maxR = k;
308  }
309  }
310  if (maxBinR < fHTCutR) return;
311 
312  // Find MAX bin XY and compare it with CUT
313  int maxBinXY = -1, maxXY = -1;
314  size = fHist.size();
315  for (unsigned int k = 0; k < size; k++){
316  if (fHist[k] > maxBinXY){
317  maxBinXY = fHist[k];
318  maxXY = k;
319  }
320  }
321  if (maxBinXY < fHTCut) return;
322 
323  CbmRichRingLight* ring1 = new CbmRichRingLight();
324 
325  // Find Preliminary Xc, Yc, R
326  float xc, yc, r;
327  float rx, ry, dr;
328  xc = (maxXY / fNofBinsX + 0.5f)* fDx + fCurMinX;
329  yc = (maxXY % fNofBinsX + 0.5f)* fDy + fCurMinY;
330  r = (maxR + 0.5f) * fDr;
331  for (int j = indmin; j < indmax + 1; j++) {
332  rx = fData[j].fHit.fX - xc;
333  ry = fData[j].fHit.fY - yc;
334 
335  dr = fabs(sqrt(rx * rx + ry * ry) - r);
336  if (dr > 0.6f) continue;
337  ring1->AddHit(fData[j].fHit, fData[j].fId);
338  }
339 
340  if (ring1->GetNofHits() < 7) {
341  delete ring1;
342  return;
343  }
344 
345  fFitCOP->DoFit(ring1);
346  float drCOPCut = fRmsCoeffCOP * sqrt(ring1->GetChi2() / ring1->GetNofHits());
347  if (drCOPCut > fMaxCutCOP) drCOPCut = fMaxCutCOP;
348 
349  xc = ring1->GetCenterX();
350  yc = ring1->GetCenterY();
351  r = ring1->GetRadius();
352 
353  delete ring1;
354 
355  CbmRichRingLight* ring2 = new CbmRichRingLight();
356  for (int j = indmin; j < indmax + 1; j++) {
357  rx = fData[j].fHit.fX - xc;
358  ry = fData[j].fHit.fY - yc;
359 
360  dr = fabs(sqrt(rx * rx + ry * ry) - r);
361  if (dr > drCOPCut) continue;
362  //fData[j+indmin].fIsUsed = true;
363  ring2->AddHit(fData[j].fHit, fData[j].fId);
364  }
365 
366  if (ring2->GetNofHits() < 7) {
367  delete ring2;
368  return;
369  }
370 
371  fFitCOP->DoFit(ring2);
372 
373  fANNSelect->DoSelect(ring2);
374  float select = ring2->GetSelectionNN();
375 
376  // Remove found hits only for good quality rings
377  if (select > fAnnCut) {
378  RemoveHitsAroundRing(indmin, indmax, ring2);
379  }
380 
381  if (select > -0.7) {
382  fFoundRings.push_back(ring2);
383  }
384  ring2=NULL;
385  delete ring2;
386 }
387 
389  int indmin,
390  int indmax,
391  CbmRichRingLight* ring)
392 {
393  float rms = sqrt(ring->GetChi2() / ring->GetNofHits());
394  float dCut = fRmsCoeffEl * rms;
395  if (dCut > fMaxCutEl) dCut = fMaxCutEl;
396 
397  for (int j = indmin; j < indmax + 1; j++) {
398  float rx = fData[j].fHit.fX - ring->GetCenterX();
399  float ry = fData[j].fHit.fY - ring->GetCenterY();
400 
401  float dr = fabs(sqrt(rx * rx + ry * ry) - ring->GetRadius());
402  if (dr < dCut) {
403  fData[j].fIsUsed = true;
404  }
405  }
406 }
407 
409 {
410  int nofRings = fFoundRings.size();
411  sort(fFoundRings.begin(), fFoundRings.end(), CbmRichRingComparatorMore());
412  set<unsigned short> usedHitsAll;
413  vector<unsigned short> goodRingIndex;
414  goodRingIndex.reserve(nofRings);
415  CbmRichRingLight* ring2;
416 
417  for (int iRing = 0; iRing < nofRings; iRing++){
418  CbmRichRingLight* ring = fFoundRings[iRing];
419  ring->SetRecFlag(-1);
420  int nofHits = ring->GetNofHits();
421  bool isGoodRingAll = true;
422  int nofUsedHitsAll = 0;
423  for(int iHit = 0; iHit < nofHits; iHit++){
424  set<unsigned short>::iterator it = usedHitsAll.find(ring->GetHitId(iHit));
425  if(it != usedHitsAll.end()){
426  nofUsedHitsAll++;
427  }
428  }
429  if ((float)nofUsedHitsAll/(float)nofHits > fUsedHitsAllCut){
430  isGoodRingAll = false;
431  }
432 
433  if (isGoodRingAll){
434  fFoundRings[iRing]->SetRecFlag(1);
435  for (int iRSet = 0; iRSet < goodRingIndex.size(); iRSet++){
436  ReAssignSharedHits(goodRingIndex[iRSet],iRing);
437  }
438  goodRingIndex.push_back(iRing);
439  for(int iHit = 0; iHit < nofHits; iHit++){
440  usedHitsAll.insert(ring->GetHitId(iHit));
441  }
442  }// isGoodRing
443  }// iRing
444 
445 // usedHits.clear();
446  usedHitsAll.clear();
447  goodRingIndex.clear();
448 }
449 
451  int ringInd1,
452  int ringInd2)
453 {
454  CbmRichRingLight* ring1 = fFoundRings[ringInd1];
455  CbmRichRingLight* ring2 = fFoundRings[ringInd2];
456  int nofHits1 = ring1->GetNofHits();
457  int nofHits2 = ring2->GetNofHits();
458 
459  for(int iHit1 = 0; iHit1 < nofHits1; iHit1++){
460  unsigned short hitInd1 = ring1->GetHitId(iHit1);
461  for(int iHit2 = 0; iHit2 < nofHits2; iHit2++){
462  unsigned short hitInd2 = ring2->GetHitId(iHit2);
463  if(hitInd1 != hitInd2) continue;
464  int hitIndData = GetHitIndex(hitInd1);
465  float hitX = fData[hitIndData].fHit.fX;
466  float hitY = fData[hitIndData].fHit.fY;
467  float rx1 = hitX - ring1->GetCenterX();
468  float ry1 = hitY - ring1->GetCenterY();
469  float dr1 = fabs(sqrt(rx1 * rx1 + ry1 * ry1) - ring1->GetRadius());
470 
471  float rx2 = hitX - ring2->GetCenterX();
472  float ry2 = hitY - ring2->GetCenterY();
473  float dr2 = fabs(sqrt(rx2 * rx2 + ry2 * ry2) - ring2->GetRadius());
474 
475  if (dr1 > dr2){
476  ring1->RemoveHit(hitInd1);
477  } else {
478  ring2->RemoveHit(hitInd2);
479  }
480  }//iHit2
481  }//iHit1
482 }
483 
485  unsigned short hitInd)
486 {
487  unsigned short size = fData.size();
488  for (unsigned short i = 0; i < size; i++){
489  if (fData[i].fId == hitInd) return i;
490  }
491  return -1;
492 }
493 
495  float x[],
496  float y[],
497  float *xc,
498  float *yc,
499  float *r)
500 {
501  register float t1, t2, t3, t4, t5, t6, t8, t9,
502  t10, t11, t14, t16, t19, t21, t41;
503 
504  t1 = x[1] * x[1];
505  t2 = x[2] * x[2];
506  t3 = y[1] * y[1];
507  t4 = y[2] * y[2];
508  t5 = t1 - t2 + t3 - t4;
509  t6 = y[0] - y[1];
510  t8 = x[0] * x[0];
511  t9 = y[0] * y[0];
512  t10 = t8 - t1 + t9 - t3;
513  t11 = y[1] - y[2];
514  t14 = x[1] - x[2];
515  t16 = x[0] - x[1];
516  t19 = 1.0f / (t14 * t6 - t16 * t11);
517 
518  *xc = 0.5e0 * (t5 * t6 - t10 * t11) * t19;
519  *yc = 0.5e0 * (t10 * t14 - t5 * t16) * t19;
520 
521  t21 = (x[0] - *xc)*(x[0] - *xc);
522  t41 = (y[0] - *yc)*(y[0] - *yc);
523  *r = sqrt(t21 + t41);
524 }