CMS 3D CMS Logo

CSCSegAlgoRU.cc
Go to the documentation of this file.
1 
8 #include "CSCSegAlgoRU.h"
9 #include "CSCSegFit.h"
15 #include <algorithm>
16 #include <cmath>
17 #include <iostream>
18 #include <string>
19 
21  : CSCSegmentAlgorithm(ps), myName("CSCSegAlgoRU"){
22  doCollisions = ps.getParameter<bool>("doCollisions");
23  enlarge = ps.getParameter<bool>("enlarge");
24  chi2_str_ = ps.getParameter<double>("chi2_str");
25  chi2Norm_2D_ = ps.getParameter<double>("chi2Norm_2D_");
26  dRMax = ps.getParameter<double>("dRMax");
27  dPhiMax = ps.getParameter<double>("dPhiMax");
28  dRIntMax = ps.getParameter<double>("dRIntMax");
29  dPhiIntMax = ps.getParameter<double>("dPhiIntMax");
30  chi2Max = ps.getParameter<double>("chi2Max");
31  wideSeg = ps.getParameter<double>("wideSeg");
32  minLayersApart = ps.getParameter<int>("minLayersApart");
33 
34  LogDebug("CSC") << myName << " has algorithm cuts set to: \n"
35  << "--------------------------------------------------------------------\n"
36  << "dRMax = " << dRMax << '\n'
37  << "dPhiMax = " << dPhiMax << '\n'
38  << "dRIntMax = " << dRIntMax << '\n'
39  << "dPhiIntMax = " << dPhiIntMax << '\n'
40  << "chi2Max = " << chi2Max << '\n'
41  << "wideSeg = " << wideSeg << '\n'
42  << "minLayersApart = " << minLayersApart << std::endl;
43 
44  //reset the thresholds for non-collision data
45  if(!doCollisions){
46  dRMax = 2.0;
47  dPhiMax = 2*dPhiMax;
48  dRIntMax = 2*dRIntMax;
49  dPhiIntMax = 2*dPhiIntMax;
50  chi2Norm_2D_ = 5*chi2Norm_2D_;
51  chi2_str_ = 100;
52  chi2Max = 2*chi2Max;
53  }
54 }
55 
56 std::vector<CSCSegment> CSCSegAlgoRU::buildSegments(const CSCChamber* aChamber, const ChamberHitContainer& urechits) const {
57  ChamberHitContainer rechits = urechits;
58  LayerIndex layerIndex(rechits.size());
59  int recHits_per_layer[6] = {0,0,0,0,0,0};
60  //skip events with high multiplicity of hits
61  if (rechits.size()>150){
62  return std::vector<CSCSegment>();
63  }
64  int iadd = 0;
65  for(unsigned int i = 0; i < rechits.size(); i++) {
66  recHits_per_layer[rechits[i]->cscDetId().layer()-1]++;//count rh per chamber
67  layerIndex[i] = rechits[i]->cscDetId().layer();
68  }
69  double z1 = aChamber->layer(1)->position().z();
70  double z6 = aChamber->layer(6)->position().z();
71  if (std::abs(z1) > std::abs(z6)){
72  reverse(layerIndex.begin(), layerIndex.end());
73  reverse(rechits.begin(), rechits.end());
74  }
75  if (rechits.size() < 2) {
76  return std::vector<CSCSegment>();
77  }
78  // We have at least 2 hits. We intend to try all possible pairs of hits to start
79  // segment building. 'All possible' means each hit lies on different layers in the chamber.
80  // after all same size segs are build we get rid of the overcrossed segments using the chi2 criteria
81  // the hits from the segs that are left are marked as used and are not added to segs in future iterations
82  // the hits from 3p segs are marked as used separately in order to try to assamble them in longer segments
83  // in case there is a second pass
84  // Choose first hit (as close to IP as possible) h1 and second hit
85  // (as far from IP as possible) h2 To do this we iterate over hits
86  // in the chamber by layer - pick two layers. Then we
87  // iterate over hits within each of these layers and pick h1 and h2
88  // these. If they are 'close enough' we build an empty
89  // segment. Then try adding hits to this segment.
90  // Initialize flags that a given hit has been allocated to a segment
91  BoolContainer used(rechits.size(), false);
92  BoolContainer used3p(rechits.size(), false);
93  // This is going to point to fits to hits, and its content will be used to create a CSCSegment
94  AlgoState aState;
95  aState.aChamber = aChamber;
96  aState.doCollisions = doCollisions;
97  aState.enlarge = enlarge;
98  aState.dRMax = dRMax;
99  aState.dPhiMax = dPhiMax;
100  aState.dRIntMax = dRIntMax;
101  aState.dPhiIntMax = dPhiIntMax;
102  aState.chi2Norm_2D_ = chi2Norm_2D_;
103  aState.chi2_str_ = chi2_str_;
104  aState.chi2Max = chi2Max;
105 
106  int scale_factor = 1;
107  if(aState.enlarge) scale_factor = 2;
108 
109  // Define buffer for segments we build
110  std::vector<CSCSegment> segments;
111  ChamberHitContainerCIt ib = rechits.begin();
112  ChamberHitContainerCIt ie = rechits.end();
113  // Possibly allow 3 passes, second widening scale factor for cuts, third for segments from displaced vertices
114  aState.windowScale = 1.; // scale factor for cuts
115  bool search_disp = false;
116  aState.strip_iadd = 1;
117  aState.chi2D_iadd = 1;
118  int npass = (wideSeg > 1.)? 3 : 2;
119  for (int ipass = 0; ipass < npass; ++ipass) {
120  if(aState.windowScale >1.){
121  iadd = 1;
122  aState.strip_iadd = 2*scale_factor;
123  aState.chi2D_iadd = 2*scale_factor;
124  if(aState.enlarge){
125  aState.chi2Max = 2*chi2Max;
126  if( rechits.size() <= 12 )iadd = 0;//allow 3 hit segments for low hit multiplicity chambers
127  }
128  }
129 
130  int used_rh = 0;
131  for (ChamberHitContainerCIt i1 = ib; i1 != ie; ++i1) {
132  if(used[i1-ib])used_rh++;
133  }
134 
135  //change the tresholds if it's time to look for displaced mu segments
136  if(aState.doCollisions && search_disp && int(rechits.size()-used_rh)>2){//check if there are enough recHits left to build a segment from displaced vertices
137  aState.doCollisions = false;
138  aState.windowScale = 1.; // scale factor for cuts
139  aState.dRMax = scale_factor*2.0;
140  aState.dPhiMax = scale_factor*2*aState.dPhiMax;
141  aState.dRIntMax = scale_factor*2*aState.dRIntMax;
142  aState.dPhiIntMax = scale_factor*2*aState.dPhiIntMax;
143  aState.chi2Norm_2D_ = scale_factor*5*aState.chi2Norm_2D_;
144  aState.chi2_str_ = scale_factor*100;
145  aState.chi2Max = scale_factor*2*aState.chi2Max;
146  }else{
147  search_disp = false;//make sure the flag is off
148  }
149 
150  for(unsigned int n_seg_min = 6u; n_seg_min > 2u + iadd; --n_seg_min){
151  BoolContainer common_used(rechits.size(),false);
152  std::array<BoolContainer, 120> common_used_it = {};
153  for (unsigned int i = 0; i < common_used_it.size(); i++) {
154  common_used_it[i] = common_used;
155  }
156  ChamberHitContainer best_proto_segment[120];
157  float min_chi[120] = {9999};
158  int common_it = 0;
159  bool first_proto_segment = true;
160  for (ChamberHitContainerCIt i1 = ib; i1 != ie; ++i1) {
161  bool segok = false;
162  //skip if rh is used and the layer tat has big rh multiplicity(>25RHs)
163  if(used[i1-ib] || recHits_per_layer[int(layerIndex[i1-ib])-1]>25 || (n_seg_min == 3 && used3p[i1-ib])) continue;
164  int layer1 = layerIndex[i1-ib];
165  const CSCRecHit2D* h1 = *i1;
166  for (ChamberHitContainerCIt i2 = ie-1; i2 != i1; --i2) {
167  if(used[i2-ib] || recHits_per_layer[int(layerIndex[i2-ib])-1]>25 || (n_seg_min == 3 && used3p[i2-ib])) continue;
168  int layer2 = layerIndex[i2-ib];
169  if((abs(layer2 - layer1) + 1) < int(n_seg_min)) break;//decrease n_seg_min
170  const CSCRecHit2D* h2 = *i2;
171  if (areHitsCloseInR(aState, h1, h2) && areHitsCloseInGlobalPhi(aState, h1, h2)) {
172  aState.proto_segment.clear();
173  if (!addHit(aState, h1, layer1))continue;
174  if (!addHit(aState, h2, layer2))continue;
175  // Can only add hits if already have a segment
176  if ( aState.sfit )tryAddingHitsToSegment(aState, rechits, used, layerIndex, i1, i2);
177  segok = isSegmentGood(aState, rechits);
178  if (segok) {
179  if(aState.proto_segment.size() > n_seg_min){
180  baseline(aState, n_seg_min);
181  updateParameters(aState);
182  }
183  if(aState.sfit->chi2() > aState.chi2Norm_2D_*aState.chi2D_iadd || aState.proto_segment.size() < n_seg_min) aState.proto_segment.clear();
184  if (!aState.proto_segment.empty()) {
185  updateParameters(aState);
186  //add same-size overcrossed protosegments to the collection
187  if(first_proto_segment){
188  flagHitsAsUsed(aState, rechits, common_used_it[0]);
189  min_chi[0] = aState.sfit->chi2();
190  best_proto_segment[0] = aState.proto_segment;
191  first_proto_segment = false;
192  }else{ //for the rest of found proto_segments
193  common_it++;
194  flagHitsAsUsed(aState, rechits, common_used_it[common_it]);
195  min_chi[common_it] = aState.sfit->chi2();
196  best_proto_segment[common_it] = aState.proto_segment;
197  ChamberHitContainerCIt hi, iu, ik;
198  int iter = common_it;
199  for(iu = ib; iu != ie; ++iu) {
200  for(hi = aState.proto_segment.begin(); hi != aState.proto_segment.end(); ++hi) {
201  if(*hi == *iu) {
202  int merge_nr = -1;
203  for(int k = 0; k < iter+1; k++){
204  if(common_used_it[k][iu-ib] == true){
205  if(merge_nr != -1){
206  //merge the k and merge_nr collections of flaged hits into the merge_nr collection and unmark the k collection hits
207  for(ik = ib; ik != ie; ++ik) {
208  if(common_used_it[k][ik-ib] == true){
209  common_used_it[merge_nr][ik-ib] = true;
210  common_used_it[k][ik-ib] = false;
211  }
212  }
213  //change best_protoseg if min_chi_2 is smaller
214  if(min_chi[k] < min_chi[merge_nr]){
215  min_chi[merge_nr] = min_chi[k];
216  best_proto_segment[merge_nr] = best_proto_segment[k];
217  best_proto_segment[k].clear();
218  min_chi[k] = 9999;
219  }
220  common_it--;
221  }
222  else{
223  merge_nr = k;
224  }
225  }//if(common_used[k][iu-ib] == true)
226  }//for k
227  }//if
228  }//for proto_seg
229  }//for rec_hits
230  }//else
231  }//proto seg not empty
232  }
233  } // h1 & h2 close
234  if (segok)
235  break;
236  } // i2
237  } // i1
238 
239 
240  //add the reconstructed segments
241  for(int j = 0;j < common_it+1; j++){
242  aState.proto_segment = best_proto_segment[j];
243  best_proto_segment[j].clear();
244  //SKIP empty proto-segments
245  if(aState.proto_segment.empty()) continue;
246  updateParameters(aState);
247  // Create an actual CSCSegment - retrieve all info from the fit
248  CSCSegment temp(aState.sfit->hits(), aState.sfit->intercept(),
249  aState.sfit->localdir(), aState.sfit->covarianceMatrix(), aState.sfit->chi2());
250  aState.sfit = nullptr;
251  segments.push_back(temp);
252  //if the segment has 3 hits flag them as used in a particular way
253  if(aState.proto_segment.size() == 3){
254  flagHitsAsUsed(aState, rechits, used3p);
255  }
256  else{
257  flagHitsAsUsed(aState, rechits, used);
258  }
259  aState.proto_segment.clear();
260  }
261  }//for n_seg_min
262 
263  if(search_disp){
264  //reset params and flags for the next chamber
265  search_disp = false;
266  aState.doCollisions = true;
267  aState.dRMax = 2.0;
268  aState.chi2_str_ = 100;
269  aState.dPhiMax = 0.5*aState.dPhiMax/scale_factor;
270  aState.dRIntMax = 0.5*aState.dRIntMax/scale_factor;
271  aState.dPhiIntMax = 0.5*aState.dPhiIntMax/scale_factor;
272  aState.chi2Norm_2D_ = 0.2*aState.chi2Norm_2D_/scale_factor;
273  aState.chi2Max = 0.5*aState.chi2Max/scale_factor;
274  }
275 
276  std::vector<CSCSegment>::iterator it =segments.begin();
277  bool good_segs = false;
278  while(it != segments.end()) {
279  if ((*it).nRecHits() > 3){
280  good_segs = true;
281  break;
282  }
283  ++it;
284  }
285  if (good_segs && aState.doCollisions) { // only change window if not enough good segments were found (bool can be changed to int if a >0 number of good segs is required)
286  search_disp = true;
287  continue;//proceed to search the segs from displaced vertices
288  }
289 
290  // Increase cut windows by factor of wideSeg only for collisions
291  if(!aState.doCollisions && !search_disp) break;
292  aState.windowScale = wideSeg;
293  } // ipass
294 
295  //get rid of enchansed 3p segments
296  std::vector<CSCSegment>::iterator it =segments.begin();
297  while(it != segments.end()) {
298  if((*it).nRecHits() == 3){
299  bool found_common = false;
300  const std::vector<CSCRecHit2D>& theseRH = (*it).specificRecHits();
301  for (ChamberHitContainerCIt i1 = ib; i1 != ie; ++i1) {
302  if(used[i1-ib] && used3p[i1-ib]){
303  const CSCRecHit2D* sh1 = *i1;
304  CSCDetId id = sh1->cscDetId();
305  int sh1layer = id.layer();
306  int RH_centerid = sh1->nStrips()/2;
307  int RH_centerStrip = sh1->channels(RH_centerid);
308  int RH_wg = sh1->hitWire();
309  std::vector<CSCRecHit2D>::const_iterator sh;
310  for(sh = theseRH.begin(); sh != theseRH.end(); ++sh){
311  CSCDetId idRH = sh->cscDetId();
312  //find segment hit coord
313  int shlayer = idRH.layer();
314  int SegRH_centerid = sh->nStrips()/2;
315  int SegRH_centerStrip = sh->channels(SegRH_centerid);
316  int SegRH_wg = sh->hitWire();
317  if(sh1layer == shlayer && SegRH_centerStrip == RH_centerStrip && SegRH_wg == RH_wg){
318  //remove the enchansed 3p segment
319  segments.erase(it,(it+1));
320  found_common = true;
321  break;
322  }
323  }//theserh
324  }
325  if(found_common) break;//current seg has already been erased
326  }//camber hits
327  if(!found_common)++it;
328  }//its a 3p seg
329  else{
330  ++it;//go to the next seg
331  }
332  }//while
333  // Give the segments to the CSCChamber
334  return segments;
335 }//build segments
336 
338  const BoolContainer& used, const LayerIndex& layerIndex,
339  const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2) const {
340  // Iterate over the layers with hits in the chamber
341  // Skip the layers containing the segment endpoints
342  // Test each hit on the other layers to see if it is near the segment
343  // If it is, see whether there is already a hit on the segment from the same layer
344  // - if so, and there are more than 2 hits on the segment, copy the segment,
345  // replace the old hit with the new hit. If the new segment chi2 is better
346  // then replace the original segment with the new one (by swap)
347  // - if not, copy the segment, add the hit. If the new chi2/dof is still satisfactory
348  // then replace the original segment with the new one (by swap)
349  ChamberHitContainerCIt ib = rechits.begin();
350  ChamberHitContainerCIt ie = rechits.end();
351  for (ChamberHitContainerCIt i = ib; i != ie; ++i) {
352  int layer = layerIndex[i-ib];
353  if (hasHitOnLayer(aState, layer) && aState.proto_segment.size() <= 2)continue;
354  if (layerIndex[i-ib] == layerIndex[i1-ib] || layerIndex[i-ib] == layerIndex[i2-ib] || used[i-ib])continue;
355 
356  const CSCRecHit2D* h = *i;
357  if (isHitNearSegment(aState, h)) {
358  // Don't consider alternate hits on layers holding the two starting points
359  if (hasHitOnLayer(aState, layer)) {
360  if (aState.proto_segment.size() <= 2)continue;
361  compareProtoSegment(aState, h, layer);
362  }
363  else{
364  increaseProtoSegment(aState, h, layer, aState.chi2D_iadd);
365  }
366  } // h & seg close
367  } // i
368 }
369 
370 bool CSCSegAlgoRU::areHitsCloseInR(const AlgoState& aState, const CSCRecHit2D* h1, const CSCRecHit2D* h2) const {
371  float maxWG_width[10] = {0, 0, 4.1, 5.69, 2.49, 5.06, 2.49, 5.06, 1.87, 5.06};
372  CSCDetId id = h1->cscDetId();
373  int iStn = id.iChamberType()-1;
374  //find maxWG_width for ME11 (tilt = 29deg)
375  int wg_num = h2->hitWire();
376  if(iStn == 0 || iStn == 1){
377  if (wg_num == 1){
378  maxWG_width[0] = 9.25;
379  maxWG_width[1] = 9.25;
380  }
381  if (wg_num > 1 && wg_num < 48){
382  maxWG_width[0] = 3.14;
383  maxWG_width[1] = 3.14;
384  }
385  if (wg_num == 48){
386  maxWG_width[0] = 10.75;
387  maxWG_width[1] = 10.75;
388  }
389  }
390  const CSCLayer* l1 = aState.aChamber->layer(h1->cscDetId().layer());
391  GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
392  const CSCLayer* l2 = aState.aChamber->layer(h2->cscDetId().layer());
393  GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
394  //find z to understand the direction
395  float h1z = gp1.z();
396  float h2z = gp2.z();
397  //switch off the IP check for non collisions case
398  if (!aState.doCollisions){
399  h1z = 1;
400  h2z = 1;
401  }
402 
403  if(aState.enlarge) {
404 
405  return (gp2.perp() > ((gp1.perp() - aState.dRMax*aState.strip_iadd*maxWG_width[iStn])*h2z)/h1z && gp2.perp() < ((gp1.perp() + aState.dRMax*aState.strip_iadd*maxWG_width[iStn])*h2z)/h1z)? true:false;
406 
407  }else{
408 
409  return (gp2.perp() > ((gp1.perp() - aState.dRMax*maxWG_width[iStn])*h2z)/h1z && gp2.perp() < ((gp1.perp() + aState.dRMax*maxWG_width[iStn])*h2z)/h1z)? true:false;
410 
411  }
412 }
413 
414 bool CSCSegAlgoRU::areHitsCloseInGlobalPhi(const AlgoState& aState, const CSCRecHit2D* h1, const CSCRecHit2D* h2) const {
415  float strip_width[10] = {0.003878509, 0.002958185, 0.002327105, 0.00152552, 0.00465421, 0.002327105, 0.00465421, 0.002327105, 0.00465421, 0.002327105};//in rad
416  const CSCLayer* l1 = aState.aChamber->layer(h1->cscDetId().layer());
417  GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
418  const CSCLayer* l2 = aState.aChamber->layer(h2->cscDetId().layer());
419  GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
420  float err_stpos_h1 = h1->errorWithinStrip();
421  float err_stpos_h2 = h2->errorWithinStrip();
422  CSCDetId id = h1->cscDetId();
423  int iStn = id.iChamberType()-1;
424  float dphi_incr = 0;
425  if(err_stpos_h1>0.25*strip_width[iStn] || err_stpos_h2>0.25*strip_width[iStn])dphi_incr = 0.5*strip_width[iStn];
426  float dphi12 = deltaPhi(gp1.barePhi(),gp2.barePhi());
427  return (fabs(dphi12) < (aState.dPhiMax*aState.strip_iadd+dphi_incr))? true:false; // +v
428 }
429 
430 bool CSCSegAlgoRU::isHitNearSegment(const AlgoState& aState, const CSCRecHit2D* h) const {
431  // Is hit near segment?
432  // Requires deltaphi and deltaR within ranges specified in parameter set.
433  // Note that to make intuitive cuts on delta(phi) one must work in
434  // phi range (-pi, +pi] not [0, 2pi)
435  float strip_width[10] = {0.003878509, 0.002958185, 0.002327105, 0.00152552, 0.00465421, 0.002327105, 0.00465421, 0.002327105, 0.00465421, 0.002327105};//in rad
436  const CSCLayer* l1 = aState.aChamber->layer((*(aState.proto_segment.begin()))->cscDetId().layer());
437  GlobalPoint gp1 = l1->toGlobal((*(aState.proto_segment.begin()))->localPosition());
438  const CSCLayer* l2 = aState.aChamber->layer((*(aState.proto_segment.begin()+1))->cscDetId().layer());
439  GlobalPoint gp2 = l2->toGlobal((*(aState.proto_segment.begin()+1))->localPosition());
440  float err_stpos_h1 = (*(aState.proto_segment.begin()))->errorWithinStrip();
441  float err_stpos_h2 = (*(aState.proto_segment.begin()+1))->errorWithinStrip();
442  const CSCLayer* l = aState.aChamber->layer(h->cscDetId().layer());
444  float err_stpos_h = h->errorWithinStrip();
445  float hphi = hp.phi(); // in (-pi, +pi]
446  if (hphi < 0.)
447  hphi += 2.*M_PI; // into range [0, 2pi)
448  float sphi = phiAtZ(aState, hp.z()); // in [0, 2*pi)
449  float phidif = sphi-hphi;
450  if (phidif < 0.)
451  phidif += 2.*M_PI; // into range [0, 2pi)
452  if (phidif > M_PI)
453  phidif -= 2.*M_PI; // into range (-pi, pi]
454  SVector6 r_glob;
455  CSCDetId id = h->cscDetId();
456  int iStn = id.iChamberType()-1;
457  float dphi_incr = 0;
458  float pos_str = 1;
459  //increase dPhi cut if the hit is on the edge of the strip
460  float stpos = (*h).positionWithinStrip();
461  bool centr_str = false;
462  if(iStn != 0 && iStn != 1){
463  if (stpos > -0.25 && stpos < 0.25) centr_str = true;
464  }
465  if(err_stpos_h1<0.25*strip_width[iStn] || err_stpos_h2<0.25*strip_width[iStn] || err_stpos_h < 0.25*strip_width[iStn]){
466  dphi_incr = 0.5*strip_width[iStn];
467  }else{
468  if(centr_str) pos_str = 1.3;
469  }
470  r_glob((*(aState.proto_segment.begin()))->cscDetId().layer()-1) = gp1.perp();
471  r_glob((*(aState.proto_segment.begin()+1))->cscDetId().layer()-1) = gp2.perp();
472  float R = hp.perp();
473  int layer = h->cscDetId().layer();
474  float r_interpolated = fit_r_phi(aState, r_glob,layer);
475  float dr = fabs(r_interpolated - R);
476  float maxWG_width[10] = {0, 0, 4.1, 5.69, 2.49, 5.06, 2.49, 5.06, 1.87, 5.06};
477  //find maxWG_width for ME11 (tilt = 29deg)
478  int wg_num = h->hitWire();
479  if(iStn == 0 || iStn == 1){
480  if (wg_num == 1){
481  maxWG_width[0] = 9.25;
482  maxWG_width[1] = 9.25;
483  }
484  if (wg_num > 1 && wg_num < 48){
485  maxWG_width[0] = 3.14;
486  maxWG_width[1] = 3.14;
487  }
488  if (wg_num == 48){
489  maxWG_width[0] = 10.75;
490  maxWG_width[1] = 10.75;
491  }
492  }
493 
494  if(aState.enlarge) {
495 
496  return (fabs(phidif) < aState.dPhiIntMax*aState.strip_iadd*pos_str+dphi_incr && fabs(dr) < aState.dRIntMax*aState.strip_iadd*maxWG_width[iStn])? true:false;
497 
498  }else{
499 
500  return (fabs(phidif) < aState.dPhiIntMax*aState.strip_iadd*pos_str+dphi_incr && fabs(dr) < aState.dRIntMax*maxWG_width[iStn])? true:false;
501 
502  }
503 }
504 
505 float CSCSegAlgoRU::phiAtZ(const AlgoState& aState, float z) const {
506  if ( !aState.sfit ) return 0.;
507  // Returns a phi in [ 0, 2*pi )
508  const CSCLayer* l1 = aState.aChamber->layer((*(aState.proto_segment.begin()))->cscDetId().layer());
509  GlobalPoint gp = l1->toGlobal(aState.sfit->intercept());
510  GlobalVector gv = l1->toGlobal(aState.sfit->localdir());
511  float x = gp.x() + (gv.x()/gv.z())*(z - gp.z());
512  float y = gp.y() + (gv.y()/gv.z())*(z - gp.z());
513  float phi = atan2(y, x);
514  if (phi < 0.f ) phi += 2. * M_PI;
515  return phi ;
516 }
517 
518 bool CSCSegAlgoRU::isSegmentGood(const AlgoState& aState, const ChamberHitContainer& rechitsInChamber) const {
519  // If the chamber has 20 hits or fewer, require at least 3 hits on segment
520  // If the chamber has >20 hits require at least 4 hits
521  // If it's the second cycle of the builder and there are <= 12 hits in chamber, require at least 3 hits on segment
522  //@@ THESE VALUES SHOULD BECOME PARAMETERS?
523  bool ok = false;
524  unsigned int iadd = ( rechitsInChamber.size() > 20)? 1 : 0;
525  if (aState.windowScale > 1.) {
526  iadd = 1;
527  if( rechitsInChamber.size() <= 12 && aState.enlarge) iadd = 0;
528  }
529  if (aState.proto_segment.size() >= 3+iadd)
530  ok = true;
531  return ok;
532 }
533 
534 void CSCSegAlgoRU::flagHitsAsUsed(const AlgoState& aState, const ChamberHitContainer& rechitsInChamber,
535  BoolContainer& used ) const {
536  // Flag hits on segment as used
537  ChamberHitContainerCIt ib = rechitsInChamber.begin();
539  for(hi = aState.proto_segment.begin(); hi != aState.proto_segment.end(); ++hi) {
540  for(iu = ib; iu != rechitsInChamber.end(); ++iu) {
541  if(*hi == *iu)
542  used[iu-ib] = true;
543  }
544  }
545 }
546 
547 bool CSCSegAlgoRU::addHit(AlgoState& aState, const CSCRecHit2D* aHit, int layer) const {
548  // Return true if hit was added successfully
549  // (and then parameters are updated).
550  // Return false if there is already a hit on the same layer, or insert failed.
551  ChamberHitContainer::const_iterator it;
552  for(it = aState.proto_segment.begin(); it != aState.proto_segment.end(); it++)
553  if (((*it)->cscDetId().layer() == layer) && (aHit != (*it)))
554  return false;
555  aState.proto_segment.push_back(aHit);
556  // make a fit
557  updateParameters(aState);
558  return true;
559 }
560 
562  // Delete input CSCSegFit, create a new one and make the fit
563  // delete sfit;
564  aState.sfit.reset(new CSCSegFit( aState.aChamber, aState.proto_segment ));
565  aState.sfit->fit();
566 }
567 
568 float CSCSegAlgoRU::fit_r_phi(const AlgoState& aState, const SVector6& points, int layer) const{
569  //find R or Phi on the given layer using the given points for the interpolation
570  float Sx = 0;
571  float Sy = 0;
572  float Sxx = 0;
573  float Sxy = 0;
574  for (int i=1;i<7;i++){
575  if (points(i-1)== 0.) continue;
576  Sy = Sy + (points(i-1));
577  Sx = Sx + i;
578  Sxx = Sxx + (i*i);
579  Sxy = Sxy + ((i)*points(i-1));
580  }
581  float delta = 2*Sxx - Sx*Sx;
582  float intercept = (Sxx*Sy - Sx*Sxy)/delta;
583  float slope = (2*Sxy - Sx*Sy)/delta;
584  return (intercept + slope*layer);
585 }
586 
587 void CSCSegAlgoRU::baseline(AlgoState& aState, int n_seg_min) const {
588  int nhits = aState.proto_segment.size();
589  //initialise vectors for strip position and error within strip
590  SVector6 sp;
591  SVector6 se;
592  unsigned int init_size = aState.proto_segment.size();
594  buffer.clear();
595  buffer.reserve(init_size);
596  while (buffer.size()< init_size){
597  ChamberHitContainer::iterator min;
598  int min_layer = 10;
599  for(ChamberHitContainer::iterator k = aState.proto_segment.begin(); k != aState.proto_segment.end(); k++){
600  const CSCRecHit2D* iRHk = *k;
601  CSCDetId idRHk = iRHk->cscDetId();
602  int kLayer = idRHk.layer();
603  if(kLayer < min_layer){
604  min_layer = kLayer;
605  min = k;
606  }
607  }
608  buffer.push_back(*min);
609  aState.proto_segment.erase(min);
610  }//while
611 
612  aState.proto_segment.clear();
613  for (ChamberHitContainer::const_iterator cand = buffer.begin(); cand != buffer.end(); cand++) {
614  aState.proto_segment.push_back(*cand);
615  }
616 
617  for(ChamberHitContainer::const_iterator iRH = aState.proto_segment.begin(); iRH != aState.proto_segment.end(); iRH++){
618  const CSCRecHit2D* iRHp = *iRH;
619  CSCDetId idRH = iRHp->cscDetId();
620  int kRing = idRH.ring();
621  int kStation = idRH.station();
622  int kLayer = idRH.layer();
623  // Find the strip containing this hit
624  int centerid = iRHp->nStrips()/2;
625  int centerStrip = iRHp->channels(centerid);
626  float stpos = (*iRHp).positionWithinStrip();
627  se(kLayer-1) = (*iRHp).errorWithinStrip();
628  // Take into account half-strip staggering of layers (ME1/1 has no staggering)
629  if (kStation == 1 && (kRing == 1 || kRing == 4)) sp(kLayer-1) = stpos + centerStrip;
630  else{
631  if (kLayer == 1 || kLayer == 3 || kLayer == 5) sp(kLayer-1) = stpos + centerStrip;
632  if (kLayer == 2 || kLayer == 4 || kLayer == 6) sp(kLayer-1) = stpos - 0.5 + centerStrip;
633  }
634  }
635  float chi2_str;
636  fitX(aState, sp, se, -1, -1, chi2_str);
637 
638  //-----------------------------------------------------
639  // Optimal point rejection method
640  //-----------------------------------------------------
641  float minSum = 1000;
642  int i1b = 0;
643  int i2b = 0;
644  int iworst = -1;
645  int bad_layer = -1;
646  ChamberHitContainer::const_iterator rh_to_be_deleted_1;
647  ChamberHitContainer::const_iterator rh_to_be_deleted_2;
648  if ( (chi2_str) > aState.chi2_str_*aState.chi2D_iadd){
649  for (ChamberHitContainer::const_iterator i1 = aState.proto_segment.begin(); i1 != aState.proto_segment.end();++i1) {
650  ++i1b;
651  const CSCRecHit2D* i1_1 = *i1;
652  CSCDetId idRH1 = i1_1->cscDetId();
653  int z1 = idRH1.layer();
654  i2b = i1b;
655  for (ChamberHitContainer::const_iterator i2 = i1+1; i2 != aState.proto_segment.end(); ++i2) {
656  ++i2b;
657  const CSCRecHit2D* i2_1 = *i2;
658  CSCDetId idRH2 = i2_1->cscDetId();
659  int z2 = idRH2.layer();
660  int irej = 0;
661  for ( ChamberHitContainer::const_iterator ir = aState.proto_segment.begin(); ir != aState.proto_segment.end(); ++ir) {
662  ++irej;
663  if (ir == i1 || ir == i2) continue;
664  float dsum = 0;
665  int hit_nr = 0;
666  const CSCRecHit2D* ir_1 = *ir;
667  CSCDetId idRH = ir_1->cscDetId();
668  int worst_layer = idRH.layer();
669  for (ChamberHitContainer::const_iterator i = aState.proto_segment.begin(); i != aState.proto_segment.end(); ++i) {
670  ++hit_nr;
671  const CSCRecHit2D* i_1 = *i;
672  if (i == i1 || i == i2 || i == ir) continue;
673  float slope = (sp(z2-1)-sp(z1-1))/(z2-z1);
674  float intersept = sp(z1-1) - slope*z1;
675  CSCDetId idRH = i_1->cscDetId();
676  int z = idRH.layer();
677  float di = fabs(sp(z-1) - intersept - slope*z);
678  dsum = dsum + di;
679  }//i
680  if (dsum < minSum){
681  minSum = dsum;
682  bad_layer = worst_layer;
683  iworst = irej;
684  rh_to_be_deleted_1 = ir;
685  }
686  }//ir
687  }//i2
688  }//i1
689  fitX(aState, sp, se, bad_layer, -1, chi2_str);
690  }//if chi2prob<1.0e-4
691 
692  //find worst from n-1 hits
693  int iworst2 = -1;
694  int bad_layer2 = -1;
695  if (iworst > -1 && (nhits-1) > n_seg_min && (chi2_str) > aState.chi2_str_*aState.chi2D_iadd){
696  iworst = -1;
697  float minSum = 1000;
698  int i1b = 0;
699  int i2b = 0;
700  for (ChamberHitContainer::const_iterator i1 = aState.proto_segment.begin(); i1 != aState.proto_segment.end();++i1) {
701  ++i1b;
702  const CSCRecHit2D* i1_1 = *i1;
703  CSCDetId idRH1 = i1_1->cscDetId();
704  int z1 = idRH1.layer();
705  i2b = i1b;
706  for ( ChamberHitContainer::const_iterator i2 = i1+1; i2 != aState.proto_segment.end(); ++i2) {
707  ++i2b;
708  const CSCRecHit2D* i2_1 = *i2;
709  CSCDetId idRH2 = i2_1->cscDetId();
710  int z2 = idRH2.layer();
711  int irej = 0;
712  for ( ChamberHitContainer::const_iterator ir = aState.proto_segment.begin(); ir != aState.proto_segment.end(); ++ir) {
713  ++irej;
714  int irej2 = 0;
715  if (ir == i1 || ir == i2 ) continue;
716  const CSCRecHit2D* ir_1 = *ir;
717  CSCDetId idRH = ir_1->cscDetId();
718  int worst_layer = idRH.layer();
719  for ( ChamberHitContainer::const_iterator ir2 = aState.proto_segment.begin(); ir2 != aState.proto_segment.end(); ++ir2) {
720  ++irej2;
721  if (ir2 == i1 || ir2 == i2 || ir2 ==ir ) continue;
722  float dsum = 0;
723  int hit_nr = 0;
724  const CSCRecHit2D* ir2_1 = *ir2;
725  CSCDetId idRH = ir2_1->cscDetId();
726  int worst_layer2 = idRH.layer();
727  for ( ChamberHitContainer::const_iterator i = aState.proto_segment.begin(); i != aState.proto_segment.end(); ++i) {
728  ++hit_nr;
729  const CSCRecHit2D* i_1 = *i;
730  if (i == i1 || i == i2 || i == ir|| i == ir2 ) continue;
731  float slope = (sp(z2-1)-sp(z1-1))/(z2-z1);
732  float intersept = sp(z1-1) - slope*z1;
733  CSCDetId idRH = i_1->cscDetId();
734  int z = idRH.layer();
735  float di = fabs(sp(z-1) - intersept - slope*z);
736  dsum = dsum + di;
737  }//i
738  if (dsum < minSum){
739  minSum = dsum;
740  iworst2 = irej2;
741  iworst = irej;
742  bad_layer = worst_layer;
743  bad_layer2 = worst_layer2;
744  rh_to_be_deleted_1 = ir;
745  rh_to_be_deleted_2 = ir2;
746  }
747  }//ir2
748  }//ir
749  }//i2
750  }//i1
751  fitX(aState, sp, se, bad_layer ,bad_layer2, chi2_str);
752  }//if prob(n-1)<e-4
753 
754  //----------------------------------
755  //erase bad_hits
756  //----------------------------------
757  if( iworst2-1 >= 0 && iworst2 <= int(aState.proto_segment.size()) ) {
758  aState.proto_segment.erase( rh_to_be_deleted_2);
759  }
760  if( iworst-1 >= 0 && iworst <= int(aState.proto_segment.size()) ){
761  aState.proto_segment.erase(rh_to_be_deleted_1);
762  }
763 }
764 
765 float CSCSegAlgoRU::fitX(const AlgoState& aState, SVector6 points, SVector6 errors, int ir, int ir2, float &chi2_str) const {
766  float S = 0;
767  float Sx = 0;
768  float Sy = 0;
769  float Sxx = 0;
770  float Sxy = 0;
771  float sigma2 = 0;
772  for (int i=1;i<7;i++){
773  if (i == ir || i == ir2 || points(i-1) == 0.) continue;
774  sigma2 = errors(i-1)*errors(i-1);
775  float i1 = i - 3.5;
776  S = S + (1/sigma2);
777  Sy = Sy + (points(i-1)/sigma2);
778  Sx = Sx + ((i1)/sigma2);
779  Sxx = Sxx + (i1*i1)/sigma2;
780  Sxy = Sxy + (((i1)*points(i-1))/sigma2);
781  }
782  float delta = S*Sxx - Sx*Sx;
783  float intercept = (Sxx*Sy - Sx*Sxy)/delta;
784  float slope = (S*Sxy - Sx*Sy)/delta;
785  float chi_str = 0;
786  chi2_str = 0;
787  // calculate chi2_str
788  for (int i=1;i<7;i++){
789  if (i == ir || i == ir2 || points(i-1) == 0.) continue;
790  chi_str = (points(i-1) - intercept - slope*(i-3.5))/(errors(i-1));
791  chi2_str = chi2_str + chi_str*chi_str;
792  }
793  return (intercept + slope*0);
794 }
795 
796 bool CSCSegAlgoRU::hasHitOnLayer(const AlgoState& aState, int layer) const {
797  // Is there is already a hit on this layer?
799  for(it = aState.proto_segment.begin(); it != aState.proto_segment.end(); it++)
800  if ((*it)->cscDetId().layer() == layer)
801  return true;
802  return false;
803 }
804 
805 bool CSCSegAlgoRU::replaceHit(AlgoState& aState, const CSCRecHit2D* h, int layer) const {
806  // replace a hit from a layer
807  ChamberHitContainer::const_iterator it;
808  for (it = aState.proto_segment.begin(); it != aState.proto_segment.end();) {
809  if ((*it)->cscDetId().layer() == layer)
810  it = aState.proto_segment.erase(it);
811  else
812  ++it;
813  }
814  return addHit(aState, h, layer);
815 }
816 
817 void CSCSegAlgoRU::compareProtoSegment(AlgoState& aState, const CSCRecHit2D* h, int layer) const {
818  // Copy the input CSCSegFit
819  std::unique_ptr<CSCSegFit> oldfit;
820  oldfit.reset(new CSCSegFit( aState.aChamber, aState.proto_segment ));
821  oldfit->fit();
822  ChamberHitContainer oldproto = aState.proto_segment;
823 
824  // May create a new fit
825  bool ok = replaceHit(aState, h, layer);
826  if ( (aState.sfit->chi2() >= oldfit->chi2() ) || !ok ) {
827  // keep original fit
828  aState.proto_segment = oldproto;
829  aState.sfit = std::move(oldfit); // reset to the original input fit
830  }
831 }
832 
833 void CSCSegAlgoRU::increaseProtoSegment(AlgoState& aState, const CSCRecHit2D* h, int layer, int chi2_factor) const {
834  // Creates a new fit
835  std::unique_ptr<CSCSegFit> oldfit;
836  ChamberHitContainer oldproto = aState.proto_segment;
837  oldfit.reset(new CSCSegFit( aState.aChamber, aState.proto_segment ));
838  oldfit->fit();
839 
840  bool ok = addHit(aState, h, layer);
841  //@@ TEST ON ndof<=0 IS JUST TO ACCEPT nhits=2 CASE??
842  if ( !ok || ( (aState.sfit->ndof() > 0) && (aState.sfit->chi2()/aState.sfit->ndof() >= aState.chi2Max)) ) {
843  // reset to original fit
844  aState.proto_segment = oldproto;
845  aState.sfit = std::move(oldfit);
846  }
847 }
#define LogDebug(id)
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
bool areHitsCloseInR(const AlgoState &aState, const CSCRecHit2D *h1, const CSCRecHit2D *h2) const
Utility functions.
CSCDetId cscDetId() const
Definition: CSCRecHit2D.h:52
T perp() const
Definition: PV3DBase.h:72
ROOT::Math::SVector< double, 6 > SVector6
Typedefs.
Definition: CSCSegAlgoRU.h:52
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
std::vector< const CSCRecHit2D * >::const_iterator ChamberHitContainerCIt
Definition: CSCSegAlgoRU.h:56
bool areHitsCloseInGlobalPhi(const AlgoState &aState, const CSCRecHit2D *h1, const CSCRecHit2D *h2) const
float phiAtZ(const AlgoState &aState, float z) const
static const double slope[3]
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
bool replaceHit(AlgoState &aState, const CSCRecHit2D *h, int layer) const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
const std::string myName
Definition: CSCSegAlgoRU.h:148
LocalPoint localPosition() const override
Definition: CSCRecHit2D.h:50
int layer() const
Definition: CSCDetId.h:61
float dPhiIntMax
Definition: CSCSegAlgoRU.h:159
ChamberHitContainer proto_segment
Definition: CSCSegAlgoRU.h:90
std::unique_ptr< CSCSegFit > sfit
Definition: CSCSegAlgoRU.h:89
void tryAddingHitsToSegment(AlgoState &aState, const ChamberHitContainer &rechitsInChamber, const BoolContainer &used, const LayerIndex &layerIndex, const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2) const
const CSCChamber * aChamber
Definition: CSCSegAlgoRU.h:85
T barePhi() const
Definition: PV3DBase.h:68
int channels(unsigned int i) const
Extracting strip channel numbers comprising the rechit - low.
Definition: CSCRecHit2D.h:55
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:48
susybsm::HSCParticleRefProd hp
Definition: classes.h:27
T z() const
Definition: PV3DBase.h:64
unsigned int nStrips() const
Definition: CSCRecHit2D.h:56
bool isSegmentGood(const AlgoState &aState, const ChamberHitContainer &rechitsInChamber) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void baseline(AlgoState &aState, int n_seg_min) const
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:39
double f[11][100]
bool isHitNearSegment(const AlgoState &aState, const CSCRecHit2D *h) const
T min(T a, T b)
Definition: MathUtil.h:58
unsigned short iChamberType() const
Definition: CSCDetId.h:107
static const std::string kLayer("layer")
void increaseProtoSegment(AlgoState &aState, const CSCRecHit2D *h, int layer, int chi2_factor) const
CSCSegAlgoRU(const edm::ParameterSet &ps)
Constructor.
Definition: CSCSegAlgoRU.cc:20
std::vector< const CSCRecHit2D * > ChamberHitContainer
Definition: CSCSegAlgoRU.h:55
std::vector< bool > BoolContainer
Definition: CSCSegAlgoRU.h:63
std::vector< int > LayerIndex
Definition: CSCSegAlgoRU.h:54
void updateParameters(AlgoState &aState) const
#define M_PI
void compareProtoSegment(AlgoState &aState, const CSCRecHit2D *h, int layer) const
int k[5][pyjets_maxn]
int ring() const
Definition: CSCDetId.h:75
float chi2Norm_2D_
Definition: CSCSegAlgoRU.h:162
bool addHit(AlgoState &aState, const CSCRecHit2D *hit, int layer) const
Utility functions.
float fit_r_phi(const AlgoState &aState, const SVector6 &points, int layer) const
short int hitWire() const
L1A.
Definition: CSCRecHit2D.h:62
std::vector< CSCSegment > buildSegments(const CSCChamber *aChamber, const ChamberHitContainer &rechits) const
Definition: CSCSegAlgoRU.cc:56
float fitX(const AlgoState &aState, SVector6 points, SVector6 errors, int ir, int ir2, float &chi2_str) const
void flagHitsAsUsed(const AlgoState &aState, const ChamberHitContainer &rechitsInChamber, BoolContainer &used) const
const int init_size
bool hasHitOnLayer(const AlgoState &aState, int layer) const
float errorWithinStrip() const
The uncertainty of the estimated position within the strip.
Definition: CSCRecHit2D.h:79
int station() const
Definition: CSCDetId.h:86
Definition: errors.py:1
T x() const
Definition: PV3DBase.h:62
def move(src, dest)
Definition: eostools.py:511
ib
Definition: cuy.py:662