CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCSegAlgoST.cc
Go to the documentation of this file.
1 
10 #include "CSCSegAlgoST.h"
11 #include "CSCSegAlgoShowering.h"
12 
14 
15 // // For clhep Matrix::solve
16 // #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
18 
21 
23 
24 #include <algorithm>
25 #include <cmath>
26 #include <iostream>
27 #include <string>
28 
29 
30 /* Constructor
31  *
32  */
33 CSCSegAlgoST::CSCSegAlgoST(const edm::ParameterSet& ps) : CSCSegmentAlgorithm(ps), myName("CSCSegAlgoST") {
34 
35  debug = ps.getUntrackedParameter<bool>("CSCDebug");
36  // minLayersApart = ps.getParameter<int>("minLayersApart");
37  // nSigmaFromSegment = ps.getParameter<double>("nSigmaFromSegment");
38  minHitsPerSegment = ps.getParameter<int>("minHitsPerSegment");
39  // muonsPerChamberMax = ps.getParameter<int>("CSCSegmentPerChamberMax");
40  // chi2Max = ps.getParameter<double>("chi2Max");
41  dXclusBoxMax = ps.getParameter<double>("dXclusBoxMax");
42  dYclusBoxMax = ps.getParameter<double>("dYclusBoxMax");
43  preClustering = ps.getParameter<bool>("preClustering");
44  preClustering_useChaining = ps.getParameter<bool>("preClusteringUseChaining");
45  Pruning = ps.getParameter<bool>("Pruning");
46  BrutePruning = ps.getParameter<bool>("BrutePruning");
47  BPMinImprovement = ps.getParameter<double>("BPMinImprovement");
48  // maxRecHitsInCluster is the maximal number of hits in a precluster that is being processed
49  // This cut is intended to remove messy events. Currently nothing is returned if there are
50  // more that maxRecHitsInCluster hits. It could be useful to return an estimate of the
51  // cluster position, which is available.
52  maxRecHitsInCluster = ps.getParameter<int>("maxRecHitsInCluster");
53  onlyBestSegment = ps.getParameter<bool>("onlyBestSegment");
54 
55  hitDropLimit4Hits = ps.getParameter<double>("hitDropLimit4Hits");
56  hitDropLimit5Hits = ps.getParameter<double>("hitDropLimit5Hits");
57  hitDropLimit6Hits = ps.getParameter<double>("hitDropLimit6Hits");
58 
59  yweightPenaltyThreshold = ps.getParameter<double>("yweightPenaltyThreshold");
60  yweightPenalty = ps.getParameter<double>("yweightPenalty");
61 
62  curvePenaltyThreshold = ps.getParameter<double>("curvePenaltyThreshold");
63  curvePenalty = ps.getParameter<double>("curvePenalty");
64 
65  useShowering = ps.getParameter<bool>("useShowering");
66  showering_ = new CSCSegAlgoShowering( ps );
67  // std::cout<<"Constructor called..."<<std::endl;
69  correctCov_ = ps.getParameter<bool>("CorrectTheErrors");
70  chi2Norm_2D_ = ps.getParameter<double>("NormChi2Cut2D");
71  chi2Norm_3D_ = ps.getParameter<double>("NormChi2Cut3D");
72  prePrun_ = ps.getParameter<bool>("prePrun");
73  prePrunLimit_ = ps.getParameter<double>("prePrunLimit");
74  //
75  condSeed1_ = ps.getParameter<double>("SeedSmall");
76  condSeed2_ = ps.getParameter<double>("SeedBig");
77  covToAnyNumber_ = ps.getParameter<bool>("ForceCovariance");
78  covToAnyNumberAll_ = ps.getParameter<bool>("ForceCovarianceAll");
79  covAnyNumber_ = ps.getParameter<double>("Covariance");
80  passCondNumber=false;
81  passCondNumber_2=false;
83  maxContrIndex=0;
84  protoNDF = 1.;
85 
86 }
87 
88 /* Destructor
89  *
90  */
92  delete showering_;
93 }
94 
95 
96 std::vector<CSCSegment> CSCSegAlgoST::run(const CSCChamber* aChamber, ChamberHitContainer rechits) {
97 
98  // Store chamber in temp memory
99  theChamber = aChamber;
100  // pre-cluster rechits and loop over all sub clusters seperately
101  std::vector<CSCSegment> segments_temp;
102  std::vector<CSCSegment> segments;
103  std::vector<ChamberHitContainer> rechits_clusters; // this is a collection of groups of rechits
104 
105  // Define yweight penalty depending on chamber. We fixed the relative ratios, but
106  // they can be scaled by parameters:
107 
108  for(int a = 0; a<5; ++a) {
109  for(int b = 0; b<5; ++b) {
110  a_yweightPenaltyThreshold[a][b] = 0.0;
111  }
112  }
113 
123 
124  if(preClustering) {
125  // run a pre-clusterer on the given rechits to split obviously separated segment seeds:
127  // it uses X,Y,Z information; there are no configurable parameters used;
128  // the X, Y, Z "cuts" are just (much) wider than the LCT readout ones
129  // (which are actually not step functions); this new code could accomodate
130  // the clusterHits one below but we leave it for security and backward
131  // comparison reasons
132  rechits_clusters = chainHits( theChamber, rechits );
133  }
134  else{
135  // it uses X,Y information + configurable parameters
136  rechits_clusters = clusterHits( theChamber, rechits );
137  }
138  // loop over the found clusters:
139  for(std::vector<ChamberHitContainer>::iterator sub_rechits = rechits_clusters.begin(); sub_rechits != rechits_clusters.end(); ++sub_rechits ) {
140  // clear the buffer for the subset of segments:
141  segments_temp.clear();
142  // build the subset of segments:
143  segments_temp = buildSegments( (*sub_rechits) );
144  // add the found subset of segments to the collection of all segments in this chamber:
145  segments.insert( segments.end(), segments_temp.begin(), segments_temp.end() );
146  }
147  // this is the place to prune:
148  if( Pruning ) {
149  segments_temp.clear(); // segments_temp needed?!?!
150  segments_temp = prune_bad_hits( theChamber, segments );
151  segments.clear(); // segments_temp needed?!?!
152  segments.swap(segments_temp); // segments_temp needed?!?!
153  }
154  if ("ME1/a" == aChamber->specs()->chamberTypeName()){
155  findDuplicates(segments);
156  }
157  return segments;
158  }
159  else {
160  segments = buildSegments(rechits);
161  if( Pruning ) {
162  segments_temp.clear(); // segments_temp needed?!?!
163  segments_temp = prune_bad_hits( theChamber, segments );
164  segments.clear(); // segments_temp needed?!?!
165  segments.swap(segments_temp); // segments_temp needed?!?!
166  }
167  if ("ME1/a" == aChamber->specs()->chamberTypeName()){
168  findDuplicates(segments);
169  }
170  return segments;
171  //return buildSegments(rechits);
172  }
173 }
174 
175 // ********************************************************************;
176 // *** This method is meant to remove clear bad hits, using as ***;
177 // *** much information from the chamber as possible (e.g. charge, ***;
178 // *** hit position, timing, etc.) ***;
179 // ********************************************************************;
180 std::vector<CSCSegment> CSCSegAlgoST::prune_bad_hits(const CSCChamber* aChamber, std::vector<CSCSegment> & segments) {
181 
182  // std::cout<<"*************************************************************"<<std::endl;
183  // std::cout<<"Called prune_bad_hits in Chamber "<< theChamber->specs()->chamberTypeName()<<std::endl;
184  // std::cout<<"*************************************************************"<<std::endl;
185 
186  std::vector<CSCSegment> segments_temp;
187  std::vector<ChamberHitContainer> rechits_clusters; // this is a collection of groups of rechits
188 
189  const float chi2ndfProbMin = 1.0e-4;
190  bool use_brute_force = BrutePruning;
191 
192  int hit_nr = 0;
193  int hit_nr_worst = -1;
194  //int hit_nr_2ndworst = -1;
195 
196  for(std::vector<CSCSegment>::iterator it=segments.begin(); it != segments.end(); ++it) {
197 
198  // do nothing for nhit <= minHitPerSegment
199  if( (*it).nRecHits() <= minHitsPerSegment ) continue;
200 
201  if( !use_brute_force ) {// find worst hit
202 
203  float chisq = (*it).chi2();
204  int nhits = (*it).nRecHits();
205  LocalPoint localPos = (*it).localPosition();
206  LocalVector segDir = (*it).localDirection();
207  const CSCChamber* cscchamber = theChamber;
208  float globZ ;
209 
210  GlobalPoint globalPosition = cscchamber->toGlobal(localPos);
211  globZ = globalPosition.z();
212 
213 
214  if( ChiSquaredProbability((double)chisq,(double)(2*nhits-4)) < chi2ndfProbMin ) {
215 
216  // find (rough) "residuals" (NOT excluding the hit from the fit - speed!) of hits on segment
217  std::vector<CSCRecHit2D> theseRecHits = (*it).specificRecHits();
218  std::vector<CSCRecHit2D>::const_iterator iRH_worst;
219  //float xdist_local = -99999.;
220 
221  float xdist_local_worst_sig = -99999.;
222  float xdist_local_2ndworst_sig = -99999.;
223  float xdist_local_sig = -99999.;
224 
225  hit_nr = 0;
226  hit_nr_worst = -1;
227  //hit_nr_2ndworst = -1;
228 
229  for ( std::vector<CSCRecHit2D>::const_iterator iRH = theseRecHits.begin(); iRH != theseRecHits.end(); ++iRH) {
230  //mark "worst" hit:
231 
232  //float z_at_target ;
233  //float radius ;
234  float loc_x_at_target ;
235  //float loc_y_at_target ;
236  //float loc_z_at_target ;
237 
238  //z_at_target = 0.;
239  loc_x_at_target = 0.;
240  //loc_y_at_target = 0.;
241  //loc_z_at_target = 0.;
242  //radius = 0.;
243 
244  // set the z target in CMS global coordinates:
245  const CSCLayer* csclayerRH = theChamber->layer((*iRH).cscDetId().layer());
246  LocalPoint localPositionRH = (*iRH).localPosition();
247  GlobalPoint globalPositionRH = csclayerRH->toGlobal(localPositionRH);
248 
249  LocalError rerrlocal = (*iRH).localPositionError();
250  float xxerr = rerrlocal.xx();
251 
252  float target_z = globalPositionRH.z(); // target z position in cm (z pos of the hit)
253 
254  if(target_z > 0.) {
255  loc_x_at_target = localPos.x() + (segDir.x()/fabs(segDir.z())*( target_z - globZ ));
256  //loc_y_at_target = localPos.y() + (segDir.y()/fabs(segDir.z())*( target_z - globZ ));
257  //loc_z_at_target = target_z;
258  }
259  else {
260  loc_x_at_target = localPos.x() + ((-1)*segDir.x()/fabs(segDir.z())*( target_z - globZ ));
261  //loc_y_at_target = localPos.y() + ((-1)*segDir.y()/fabs(segDir.z())*( target_z - globZ ));
262  //loc_z_at_target = target_z;
263  }
264  // have to transform the segments coordinates back to the local frame... how?!!!!!!!!!!!!
265 
266  //xdist_local = fabs(localPositionRH.x() - loc_x_at_target);
267  xdist_local_sig = fabs((localPositionRH.x() -loc_x_at_target)/(xxerr));
268 
269  if( xdist_local_sig > xdist_local_worst_sig ) {
270  xdist_local_2ndworst_sig = xdist_local_worst_sig;
271  xdist_local_worst_sig = xdist_local_sig;
272  iRH_worst = iRH;
273  //hit_nr_2ndworst = hit_nr_worst;
274  hit_nr_worst = hit_nr;
275  }
276  else if(xdist_local_sig > xdist_local_2ndworst_sig) {
277  xdist_local_2ndworst_sig = xdist_local_sig;
278  //hit_nr_2ndworst = hit_nr;
279  }
280  ++hit_nr;
281  }
282 
283  // reset worst hit number if certain criteria apply.
284  // Criteria: 2nd worst hit must be at least a factor of
285  // 1.5 better than the worst in terms of sigma:
286  if ( xdist_local_worst_sig / xdist_local_2ndworst_sig < 1.5 ) {
287  hit_nr_worst = -1;
288  //hit_nr_2ndworst = -1;
289  }
290  }
291  }
292 
293  // if worst hit was found, refit without worst hit and select if considerably better than original fit.
294  // Can also use brute force: refit all n-1 hit segments and choose one over the n hit if considerably "better"
295 
296  std::vector< CSCRecHit2D > buffer;
297  std::vector< std::vector< CSCRecHit2D > > reduced_segments;
298  std::vector< CSCRecHit2D > theseRecHits = (*it).specificRecHits();
299  float best_red_seg_prob = 0.0;
300  // usefor chi2 1 diff float best_red_seg_prob = 99999.;
301  buffer.clear();
302 
303  if( ChiSquaredProbability((double)(*it).chi2(),(double)((2*(*it).nRecHits())-4)) < chi2ndfProbMin ) {
304 
305  buffer = theseRecHits;
306 
307  // Dirty switch: here one can select to refit all possible subsets or just the one without the
308  // tagged worst hit:
309  if( use_brute_force ) { // Brute force method: loop over all possible segments:
310  for(size_t bi = 0; bi < buffer.size(); ++bi) {
311  reduced_segments.push_back(buffer);
312  reduced_segments[bi].erase(reduced_segments[bi].begin()+(bi),reduced_segments[bi].begin()+(bi+1));
313  }
314  }
315  else { // More elegant but still biased: erase only worst hit
316  // Comment: There is not a very strong correlation of the worst hit with the one that one should remove...
317  if( hit_nr_worst >= 0 && hit_nr_worst <= int(buffer.size()) ) {
318  // fill segment in buffer, delete worst hit
319  buffer.erase(buffer.begin()+(hit_nr_worst),buffer.begin()+(hit_nr_worst+1));
320  reduced_segments.push_back(buffer);
321  }
322  else {
323  // only fill segment in array, do not delete anything
324  reduced_segments.push_back(buffer);
325  }
326  }
327  }
328 
329  // Loop over the subsegments and fit (only one segment if "use_brute_force" is false):
330  for(size_t iSegment=0; iSegment<reduced_segments.size(); ++iSegment) {
331  // loop over hits on given segment and push pointers to hits into protosegment
332  protoSegment.clear();
333  for(size_t m = 0; m<reduced_segments[iSegment].size(); ++m ) {
334  protoSegment.push_back(&reduced_segments[iSegment][m]);
335  }
336  passCondNumber=false;
337  passCondNumber_2 = false;
339  doSlopesAndChi2();
340  // Attempt to handle numerical instability of the fit;
341  // The same as in the build method;
342  // Preprune is not applied;
343  if(correctCov_){
345  passCondNumber = true;
346  doSlopesAndChi2();
347  }
349  passCondNumber_2=true;
350  doSlopesAndChi2();
351  }
352  }
354  // calculate error matrix
355  AlgebraicSymMatrix protoErrors = calculateError();
356  // but reorder components to match what's required by TrackingRecHit interface
357  // i.e. slopes first, then positions
358  flipErrors( protoErrors );
359  //
361 
362  // replace n hit segment with n-1 hit segment, if segment probability is BPMinImprovement better:
363  if( ( ChiSquaredProbability((double)(*it).chi2(),(double)((2*(*it).nRecHits())-4))
364  <
365  (1./BPMinImprovement)*(ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4))) ) // was (1.e-3) 081202
366 
367  &&
368  ( (ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4)))
369  > best_red_seg_prob
370  )
371  &&
372  ( (ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4))) > 1e-10 )
373  ) {
374  best_red_seg_prob = ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4));
375  // The alternative n-1 segment is much cleaner. If this segment
376  // has >= minHitsPerSegment hits exchange current n hit segment (*it)
377  // with better n-1 hit segment:
378  if( temp.nRecHits() >= minHitsPerSegment ) {
379  (*it) = temp;
380  }
381  }
382  }
383  }
384 
385  return segments;
386 
387 }
388 
389 
390 // ********************************************************************;
391 std::vector< std::vector<const CSCRecHit2D*> > CSCSegAlgoST::clusterHits(const CSCChamber* aChamber, ChamberHitContainer & rechits) {
392  theChamber = aChamber;
393 
394  std::vector<ChamberHitContainer> rechits_clusters; // this is a collection of groups of rechits
395  // const float dXclus_box_cut = 4.; // seems to work reasonably 070116
396  // const float dYclus_box_cut = 8.; // seems to work reasonably 070116
397 
398  //float dXclus = 0.0;
399  //float dYclus = 0.0;
400  float dXclus_box = 0.0;
401  float dYclus_box = 0.0;
402 
403  std::vector<const CSCRecHit2D*> temp;
404 
405  std::vector< ChamberHitContainer > seeds;
406 
407  std::vector<float> running_meanX;
408  std::vector<float> running_meanY;
409 
410  std::vector<float> seed_minX;
411  std::vector<float> seed_maxX;
412  std::vector<float> seed_minY;
413  std::vector<float> seed_maxY;
414 
415  //std::cout<<"*************************************************************"<<std::endl;
416  //std::cout<<"Called clusterHits in Chamber "<< theChamber->specs()->chamberTypeName()<<std::endl;
417  //std::cout<<"*************************************************************"<<std::endl;
418 
419  // split rechits into subvectors and return vector of vectors:
420  // Loop over rechits
421  // Create one seed per hit
422  for(unsigned int i = 0; i < rechits.size(); ++i) {
423 
424  temp.clear();
425 
426  temp.push_back(rechits[i]);
427 
428  seeds.push_back(temp);
429 
430  // First added hit in seed defines the mean to which the next hit is compared
431  // for this seed.
432 
433  running_meanX.push_back( rechits[i]->localPosition().x() );
434  running_meanY.push_back( rechits[i]->localPosition().y() );
435 
436  // set min/max X and Y for box containing the hits in the precluster:
437  seed_minX.push_back( rechits[i]->localPosition().x() );
438  seed_maxX.push_back( rechits[i]->localPosition().x() );
439  seed_minY.push_back( rechits[i]->localPosition().y() );
440  seed_maxY.push_back( rechits[i]->localPosition().y() );
441  }
442 
443  // merge clusters that are too close
444  // measure distance between final "running mean"
445  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
446 
447  for(size_t MMM = NNN+1; MMM < seeds.size(); ++MMM) {
448  if(running_meanX[MMM] == 999999. || running_meanX[NNN] == 999999. ) {
449  LogDebug("CSCSegment|CSC") << "CSCSegmentST::clusterHits: Warning: Skipping used seeds, this should happen - inform developers!";
450  // std::cout<<"We should never see this line now!!!"<<std::endl;
451  continue; //skip seeds that have been used
452  }
453 
454  // calculate cut criteria for simple running mean distance cut:
455  //dXclus = fabs(running_meanX[NNN] - running_meanX[MMM]);
456  //dYclus = fabs(running_meanY[NNN] - running_meanY[MMM]);
457 
458  // calculate minmal distance between precluster boxes containing the hits:
459  if ( running_meanX[NNN] > running_meanX[MMM] ) dXclus_box = seed_minX[NNN] - seed_maxX[MMM];
460  else dXclus_box = seed_minX[MMM] - seed_maxX[NNN];
461  if ( running_meanY[NNN] > running_meanY[MMM] ) dYclus_box = seed_minY[NNN] - seed_maxY[MMM];
462  else dYclus_box = seed_minY[MMM] - seed_maxY[NNN];
463 
464 
465  if( dXclus_box < dXclusBoxMax && dYclus_box < dYclusBoxMax ) {
466  // merge clusters!
467  // merge by adding seed NNN to seed MMM and erasing seed NNN
468 
469  // calculate running mean for the merged seed:
470  running_meanX[MMM] = (running_meanX[NNN]*seeds[NNN].size() + running_meanX[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
471  running_meanY[MMM] = (running_meanY[NNN]*seeds[NNN].size() + running_meanY[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
472 
473  // update min/max X and Y for box containing the hits in the merged cluster:
474  if ( seed_minX[NNN] <= seed_minX[MMM] ) seed_minX[MMM] = seed_minX[NNN];
475  if ( seed_maxX[NNN] > seed_maxX[MMM] ) seed_maxX[MMM] = seed_maxX[NNN];
476  if ( seed_minY[NNN] <= seed_minY[MMM] ) seed_minY[MMM] = seed_minY[NNN];
477  if ( seed_maxY[NNN] > seed_maxY[MMM] ) seed_maxY[MMM] = seed_maxY[NNN];
478 
479  // add seed NNN to MMM (lower to larger number)
480  seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
481 
482  // mark seed NNN as used (at the moment just set running mean to 999999.)
483  running_meanX[NNN] = 999999.;
484  running_meanY[NNN] = 999999.;
485  // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to
486  // next seed (NNN+1)
487  break;
488  }
489 
490  }
491  }
492 
493  // hand over the final seeds to the output
494  // would be more elegant if we could do the above step with
495  // erasing the merged ones, rather than the
496  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
497  if(running_meanX[NNN] == 999999.) continue; //skip seeds that have been marked as used up in merging
498  rechits_clusters.push_back(seeds[NNN]);
499  }
500 
501  //***************************************************************
502 
503  return rechits_clusters;
504 }
505 
506 
507 std::vector< std::vector<const CSCRecHit2D*> > CSCSegAlgoST::chainHits(const CSCChamber* aChamber, ChamberHitContainer & rechits) {
508 
509  std::vector<ChamberHitContainer> rechits_chains; // this is a collection of groups of rechits
510 
511 
512  std::vector<const CSCRecHit2D*> temp;
513 
514  std::vector< ChamberHitContainer > seeds;
515 
516  std::vector <bool> usedCluster;
517 
518  // split rechits into subvectors and return vector of vectors:
519  // Loop over rechits
520  // Create one seed per hit
521  //std::cout<<" rechits.size() = "<<rechits.size()<<std::endl;
522  for(unsigned int i = 0; i < rechits.size(); ++i) {
523  temp.clear();
524  temp.push_back(rechits[i]);
525  seeds.push_back(temp);
526  usedCluster.push_back(false);
527  }
528  bool isME11a = false;
529  if ("ME1/a" == aChamber->specs()->chamberTypeName()){
530  isME11a = true;
531  }
532  // merge chains that are too close ("touch" each other)
533  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
534  for(size_t MMM = NNN+1; MMM < seeds.size(); ++MMM) {
535  if(usedCluster[MMM] || usedCluster[NNN]){
536  continue;
537  }
538  // all is in the way we define "good";
539  // try not to "cluster" the hits but to "chain" them;
540  // it does the clustering but also does a better job
541  // for inclined tracks (not clustering them togheter;
542  // crossed tracks would be still clustered together)
543  // 22.12.09: In fact it is not much more different
544  // than the "clustering", we just introduce another
545  // variable in the game - Z. And it make sense
546  // to re-introduce Y (or actually wire group mumber)
547  // in a similar way as for the strip number - see
548  // the code below.
549  bool goodToMerge = isGoodToMerge(isME11a, seeds[NNN], seeds[MMM]);
550  if(goodToMerge){
551  // merge chains!
552  // merge by adding seed NNN to seed MMM and erasing seed NNN
553 
554  // add seed NNN to MMM (lower to larger number)
555  seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
556 
557  // mark seed NNN as used
558  usedCluster[NNN] = true;
559  // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to
560  // next seed (NNN+1)
561  break;
562  }
563 
564  }
565  }
566 
567  // hand over the final seeds to the output
568  // would be more elegant if we could do the above step with
569  // erasing the merged ones, rather than the
570 
571  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
572  if(usedCluster[NNN]) continue; //skip seeds that have been marked as used up in merging
573  rechits_chains.push_back(seeds[NNN]);
574  }
575 
576  //***************************************************************
577 
578  return rechits_chains;
579 }
580 
581 bool CSCSegAlgoST::isGoodToMerge(bool isME11a, ChamberHitContainer & newChain, ChamberHitContainer & oldChain) {
582  for(size_t iRH_new = 0;iRH_new<newChain.size();++iRH_new){
583  int layer_new = newChain[iRH_new]->cscDetId().layer()-1;
584  int middleStrip_new = newChain[iRH_new]->nStrips()/2;
585  int centralStrip_new = newChain[iRH_new]->channels(middleStrip_new);
586  int centralWire_new = newChain[iRH_new]->hitWire();
587  bool layerRequirementOK = false;
588  bool stripRequirementOK = false;
589  bool wireRequirementOK = false;
590  bool goodToMerge = false;
591  for(size_t iRH_old = 0;iRH_old<oldChain.size();++iRH_old){
592  int layer_old = oldChain[iRH_old]->cscDetId().layer()-1;
593  int middleStrip_old = oldChain[iRH_old]->nStrips()/2;
594  int centralStrip_old = oldChain[iRH_old]->channels(middleStrip_old);
595  int centralWire_old = oldChain[iRH_old]->hitWire();
596 
597  // to be chained, two hits need to be in neighbouring layers...
598  // or better allow few missing layers (upto 3 to avoid inefficiencies);
599  // however we'll not make an angle correction because it
600  // worsen the situation in some of the "regular" cases
601  // (not making the correction means that the conditions for
602  // forming a cluster are different if we have missing layers -
603  // this could affect events at the boundaries )
604  if(layer_new==layer_old+1 || layer_new==layer_old-1 ||
605  layer_new==layer_old+2 || layer_new==layer_old-2 ||
606  layer_new==layer_old+3 || layer_new==layer_old-3 ||
607  layer_new==layer_old+4 || layer_new==layer_old-4 ){
608  layerRequirementOK = true;
609  }
610  int allStrips = 48;
611  //to be chained, two hits need to be "close" in strip number (can do it in phi
612  // but it doesn't really matter); let "close" means upto 2 strips (3?) -
613  // this is more compared to what CLCT readout patterns allow
614  if(centralStrip_new==centralStrip_old ||
615  centralStrip_new==centralStrip_old+1 || centralStrip_new==centralStrip_old-1 ||
616  centralStrip_new==centralStrip_old+2 || centralStrip_new==centralStrip_old-2){
617  stripRequirementOK = true;
618  }
619  // same for wires (and ALCT patterns)
620  if(centralWire_new==centralWire_old ||
621  centralWire_new==centralWire_old+1 || centralWire_new==centralWire_old-1 ||
622  centralWire_new==centralWire_old+2 || centralWire_new==centralWire_old-2){
623  wireRequirementOK = true;
624  }
625 
626  if(isME11a){
627  if(centralStrip_new==centralStrip_old+1-allStrips || centralStrip_new==centralStrip_old-1-allStrips ||
628  centralStrip_new==centralStrip_old+2-allStrips || centralStrip_new==centralStrip_old-2-allStrips ||
629  centralStrip_new==centralStrip_old+1+allStrips || centralStrip_new==centralStrip_old-1+allStrips ||
630  centralStrip_new==centralStrip_old+2+allStrips || centralStrip_new==centralStrip_old-2+allStrips){
631  stripRequirementOK = true;
632  }
633  }
634  if(layerRequirementOK && stripRequirementOK && wireRequirementOK){
635  goodToMerge = true;
636  return goodToMerge;
637  }
638  }
639  }
640  return false;
641 }
642 
643 
644 
645 
646 double CSCSegAlgoST::theWeight(double coordinate_1, double coordinate_2, double coordinate_3, float layer_1, float layer_2, float layer_3) {
647  double sub_weight = 0;
648  sub_weight = fabs(
649  ( (coordinate_2 - coordinate_3) / (layer_2 - layer_3) ) -
650  ( (coordinate_1 - coordinate_2) / (layer_1 - layer_2) )
651  );
652  return sub_weight;
653 }
654 
655 /*
656  * This algorithm is based on the Minimum Spanning Tree (ST) approach
657  * for building endcap muon track segments out of the rechit's in a CSCChamber.
658  */
660 
661  // Clear buffer for segment vector
662  std::vector<CSCSegment> segmentInChamber;
663  segmentInChamber.clear(); // list of final segments
664 
665  // CSC Ring;
666  unsigned int thering = 999;
667  unsigned int thestation = 999;
668  //unsigned int thecham = 999;
669 
670  std::vector<int> hits_onLayerNumber(6);
671 
672  unsigned int UpperLimit = maxRecHitsInCluster;
673  if (int(rechits.size()) < minHitsPerSegment) return segmentInChamber;
674 
675  for(int iarray = 0; iarray <6; ++iarray) { // magic number 6: number of layers in CSC chamber - not gonna change :)
676  PAhits_onLayer[iarray].clear();
677  hits_onLayerNumber[iarray] = 0;
678  }
679 
680  chosen_Psegments.clear();
681  chosen_weight_A.clear();
682 
683  Psegments.clear();
684  Psegments_noLx.clear();
685  Psegments_noL1.clear();
686  Psegments_noL2.clear();
687  Psegments_noL3.clear();
688  Psegments_noL4.clear();
689  Psegments_noL5.clear();
690  Psegments_noL6.clear();
691 
692  Psegments_hits.clear();
693 
694  weight_A.clear();
695  weight_noLx_A.clear();
696  weight_noL1_A.clear();
697  weight_noL2_A.clear();
698  weight_noL3_A.clear();
699  weight_noL4_A.clear();
700  weight_noL5_A.clear();
701  weight_noL6_A.clear();
702 
703  weight_B.clear();
704  weight_noL1_B.clear();
705  weight_noL2_B.clear();
706  weight_noL3_B.clear();
707  weight_noL4_B.clear();
708  weight_noL5_B.clear();
709  weight_noL6_B.clear();
710 
711  curv_A.clear();
712  curv_noL1_A.clear();
713  curv_noL2_A.clear();
714  curv_noL3_A.clear();
715  curv_noL4_A.clear();
716  curv_noL5_A.clear();
717  curv_noL6_A.clear();
718 
719  // definition of middle layer for n-hit segment
720  int midlayer_pointer[6] = {0,0,2,3,3,4};
721 
722  // int n_layers_missed_tot = 0;
723  int n_layers_occupied_tot = 0;
724  int n_layers_processed = 0;
725 
726  float min_weight_A = 99999.9;
727  float min_weight_noLx_A = 99999.9;
728 
729  //float best_weight_B = -1.;
730  //float best_weight_noLx_B = -1.;
731 
732  //float best_curv_A = -1.;
733  //float best_curv_noLx_A = -1.;
734 
735  int best_pseg = -1;
736  int best_noLx_pseg = -1;
737  int best_Layer_noLx = -1;
738 
739  //************************************************************************;
740  //*** Start segment building *****************************************;
741  //************************************************************************;
742 
743  // Determine how many layers with hits we have
744  // Fill all hits into the layer hit container:
745 
746  // Have 2 standard arrays: one giving the number of hits per layer.
747  // The other the corresponding hits.
748 
749  // Loop all available hits, count hits per layer and fill the hits into array by layer
750  for(size_t M = 0; M < rechits.size(); ++M) {
751  // add hits to array per layer and count hits per layer:
752  hits_onLayerNumber[ rechits[M]->cscDetId().layer()-1 ] += 1;
753  if(hits_onLayerNumber[ rechits[M]->cscDetId().layer()-1 ] == 1 ) n_layers_occupied_tot += 1;
754  // add hits to vector in array
755  PAhits_onLayer[rechits[M]->cscDetId().layer()-1] .push_back(rechits[M]);
756  }
757 
758  // We have now counted the hits per layer and filled pointers to the hits into an array
759 
760  int tothits = 0;
761  int maxhits = 0;
762  int nexthits = 0;
763  int maxlayer = -1;
764  int nextlayer = -1;
765 
766  for(size_t i = 0; i< hits_onLayerNumber.size(); ++i){
767  //std::cout<<"We have "<<hits_onLayerNumber[i]<<" hits on layer "<<i+1<<std::endl;
768  tothits += hits_onLayerNumber[i];
769  if (hits_onLayerNumber[i] > maxhits) {
770  nextlayer = maxlayer;
771  nexthits = maxhits;
772  maxlayer = i;
773  maxhits = hits_onLayerNumber[i];
774  }
775  else if (hits_onLayerNumber[i] > nexthits) {
776  nextlayer = i;
777  nexthits = hits_onLayerNumber[i];
778  }
779  }
780 
781 
782  if (tothits > (int)UpperLimit) {
783  if (n_layers_occupied_tot > 4) {
784  tothits = tothits - hits_onLayerNumber[maxlayer];
785  n_layers_occupied_tot = n_layers_occupied_tot - 1;
786  PAhits_onLayer[maxlayer].clear();
787  hits_onLayerNumber[maxlayer] = 0;
788  }
789  }
790 
791  if (tothits > (int)UpperLimit) {
792  if (n_layers_occupied_tot > 4) {
793  tothits = tothits - hits_onLayerNumber[nextlayer];
794  n_layers_occupied_tot = n_layers_occupied_tot - 1;
795  PAhits_onLayer[nextlayer].clear();
796  hits_onLayerNumber[nextlayer] = 0;
797  }
798  }
799 
800  if (tothits > (int)UpperLimit){
801 
802  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
803  // Showering muon - returns nothing if chi2 == -1 (see comment in SegAlgoShowering)
804  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
805  if (useShowering) {
806  CSCSegment segShower = showering_->showerSeg(theChamber, rechits);
807 
808  // Make sure have at least 3 hits...
809  if ( segShower.nRecHits() < 3 ) return segmentInChamber;
810  if ( segShower.chi2() == -1 ) return segmentInChamber;
811 
812  segmentInChamber.push_back(segShower);
813  return segmentInChamber;
814 
815  } else{
816  LogDebug("CSCSegment|CSC") <<"Number of rechits in the cluster/chamber > "<< UpperLimit<<
817  " ... Segment finding in the cluster/chamber canceled!";
818  // std::cout<<"Number of rechits in the cluster/chamber > "<< UpperLimit<<
819  // " ... Segment finding in the cluster/chamber canceled! "<<std::endl;
820  return segmentInChamber;
821  }
822  }
823 
824  // Find out which station, ring and chamber we are in
825  // Used to choose station/ring dependant y-weight cuts
826 
827  if( rechits.size() > 0 ) {
828  thering = rechits[0]->cscDetId().ring();
829  thestation = rechits[0]->cscDetId().station();
830  //thecham = rechits[0]->cscDetId().chamber();
831  }
832 
833  // std::cout<<"We are in Station/ring/chamber: "<<thestation <<" "<< thering<<" "<< thecham<<std::endl;
834 
835  // Cut-off parameter - don't reconstruct segments with less than X hits
836  if( n_layers_occupied_tot < minHitsPerSegment ) {
837  return segmentInChamber;
838  }
839 
840  // Start building all possible hit combinations:
841 
842  // loop over the six chamber layers and form segment candidates from the available hits:
843 
844  for(int layer = 0; layer < 6; ++layer) {
845 
846  // *****************************************************************
847  // *** Set missed layer counter here (not currently implemented) ***
848  // *****************************************************************
849  // if( PAhits_onLayer[layer].size() == 0 ) {
850  // n_layers_missed_tot += 1;
851  // }
852 
853  if( PAhits_onLayer[layer].size() > 0 ) {
854  n_layers_processed += 1;
855  }
856 
857  // Save the size of the protosegment before hits were added on the current layer
858  int orig_number_of_psegs = Psegments.size();
859  int orig_number_of_noL1_psegs = Psegments_noL1.size();
860  int orig_number_of_noL2_psegs = Psegments_noL2.size();
861  int orig_number_of_noL3_psegs = Psegments_noL3.size();
862  int orig_number_of_noL4_psegs = Psegments_noL4.size();
863  int orig_number_of_noL5_psegs = Psegments_noL5.size();
864  int orig_number_of_noL6_psegs = Psegments_noL6.size();
865 
866  // loop over the hits on the layer and initiate protosegments or add hits to protosegments
867  for(int hit = 0; hit < int(PAhits_onLayer[layer].size()); ++hit) { // loop all hits on the Layer number "layer"
868 
869  // create protosegments from all hits on the first layer with hits
870  if( orig_number_of_psegs == 0 ) { // would be faster to turn this around - ask for "orig_number_of_psegs != 0"
871 
872  Psegments_hits.push_back(PAhits_onLayer[layer][hit]);
873 
874  Psegments.push_back(Psegments_hits);
875  Psegments_noL6.push_back(Psegments_hits);
876  Psegments_noL5.push_back(Psegments_hits);
877  Psegments_noL4.push_back(Psegments_hits);
878  Psegments_noL3.push_back(Psegments_hits);
879  Psegments_noL2.push_back(Psegments_hits);
880 
881  // Initialize weights corresponding to this segment for first hit (with 0)
882 
883  curv_A.push_back(0.0);
884  curv_noL6_A.push_back(0.0);
885  curv_noL5_A.push_back(0.0);
886  curv_noL4_A.push_back(0.0);
887  curv_noL3_A.push_back(0.0);
888  curv_noL2_A.push_back(0.0);
889 
890  weight_A.push_back(0.0);
891  weight_noL6_A.push_back(0.0);
892  weight_noL5_A.push_back(0.0);
893  weight_noL4_A.push_back(0.0);
894  weight_noL3_A.push_back(0.0);
895  weight_noL2_A.push_back(0.0);
896 
897  weight_B.push_back(0.0);
898  weight_noL6_B.push_back(0.0);
899  weight_noL5_B.push_back(0.0);
900  weight_noL4_B.push_back(0.0);
901  weight_noL3_B.push_back(0.0);
902  weight_noL2_B.push_back(0.0);
903 
904  // reset array for next hit on next layer
905  Psegments_hits .clear();
906  }
907  else {
908  if( orig_number_of_noL1_psegs == 0 ) {
909 
910  Psegments_hits.push_back(PAhits_onLayer[layer][hit]);
911 
912  Psegments_noL1.push_back(Psegments_hits);
913 
914  // Initialize weight corresponding to this segment for first hit (with 0)
915 
916  curv_noL1_A.push_back(0.0);
917 
918  weight_noL1_A.push_back(0.0);
919 
920  weight_noL1_B.push_back(0.0);
921 
922  // reset array for next hit on next layer
923  Psegments_hits .clear();
924 
925  }
926 
927  // loop over the protosegments and create a new protosegments for each hit-1 on this layer
928 
929  for( int pseg = 0; pseg < orig_number_of_psegs; ++pseg ) {
930 
931  int pseg_pos = (pseg)+((hit)*orig_number_of_psegs);
932  int pseg_noL1_pos = (pseg)+((hit)*orig_number_of_noL1_psegs);
933  int pseg_noL2_pos = (pseg)+((hit)*orig_number_of_noL2_psegs);
934  int pseg_noL3_pos = (pseg)+((hit)*orig_number_of_noL3_psegs);
935  int pseg_noL4_pos = (pseg)+((hit)*orig_number_of_noL4_psegs);
936  int pseg_noL5_pos = (pseg)+((hit)*orig_number_of_noL5_psegs);
937  int pseg_noL6_pos = (pseg)+((hit)*orig_number_of_noL6_psegs);
938 
939  // - Loop all psegs.
940  // - If not last hit, clone existing protosegments (PAhits_onLayer[layer].size()-1) times
941  // - then add the new hits
942 
943  if( ! (hit == int(PAhits_onLayer[layer].size()-1)) ) { // not the last hit - prepare (copy) new protosegments for the following hits
944  // clone psegs (to add next hits or last hit on layer):
945 
946  Psegments.push_back(Psegments[ pseg_pos ]);
947  if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) Psegments_noL1.push_back(Psegments_noL1[ pseg_noL1_pos ]);
948  if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) Psegments_noL2.push_back(Psegments_noL2[ pseg_noL2_pos ]);
949  if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) Psegments_noL3.push_back(Psegments_noL3[ pseg_noL3_pos ]);
950  if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) Psegments_noL4.push_back(Psegments_noL4[ pseg_noL4_pos ]);
951  if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) Psegments_noL5.push_back(Psegments_noL5[ pseg_noL5_pos ]);
952  if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) Psegments_noL6.push_back(Psegments_noL6[ pseg_noL6_pos ]);
953  // clone weight corresponding to this segment too
954  weight_A.push_back(weight_A[ pseg_pos ]);
955  if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) weight_noL1_A.push_back(weight_noL1_A[ pseg_noL1_pos ]);
956  if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) weight_noL2_A.push_back(weight_noL2_A[ pseg_noL2_pos ]);
957  if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) weight_noL3_A.push_back(weight_noL3_A[ pseg_noL3_pos ]);
958  if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) weight_noL4_A.push_back(weight_noL4_A[ pseg_noL4_pos ]);
959  if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) weight_noL5_A.push_back(weight_noL5_A[ pseg_noL5_pos ]);
960  if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) weight_noL6_A.push_back(weight_noL6_A[ pseg_noL6_pos ]);
961  // clone curvature variable corresponding to this segment too
962  curv_A.push_back(curv_A[ pseg_pos ]);
963  if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) curv_noL1_A.push_back(curv_noL1_A[ pseg_noL1_pos ]);
964  if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) curv_noL2_A.push_back(curv_noL2_A[ pseg_noL2_pos ]);
965  if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) curv_noL3_A.push_back(curv_noL3_A[ pseg_noL3_pos ]);
966  if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) curv_noL4_A.push_back(curv_noL4_A[ pseg_noL4_pos ]);
967  if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) curv_noL5_A.push_back(curv_noL5_A[ pseg_noL5_pos ]);
968  if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) curv_noL6_A.push_back(curv_noL6_A[ pseg_noL6_pos ]);
969  // clone "y"-weight corresponding to this segment too
970  weight_B.push_back(weight_B[ pseg_pos ]);
971  if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) weight_noL1_B.push_back(weight_noL1_B[ pseg_noL1_pos ]);
972  if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) weight_noL2_B.push_back(weight_noL2_B[ pseg_noL2_pos ]);
973  if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) weight_noL3_B.push_back(weight_noL3_B[ pseg_noL3_pos ]);
974  if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) weight_noL4_B.push_back(weight_noL4_B[ pseg_noL4_pos ]);
975  if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) weight_noL5_B.push_back(weight_noL5_B[ pseg_noL5_pos ]);
976  if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) weight_noL6_B.push_back(weight_noL6_B[ pseg_noL6_pos ]);
977  }
978  // add hits to original pseg:
979  Psegments[ pseg_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
980  if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) Psegments_noL1[ pseg_noL1_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
981  if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) Psegments_noL2[ pseg_noL2_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
982  if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) Psegments_noL3[ pseg_noL3_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
983  if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) Psegments_noL4[ pseg_noL4_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
984  if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) Psegments_noL5[ pseg_noL5_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
985  if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) Psegments_noL6[ pseg_noL6_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
986 
987  // calculate/update the weight (only for >2 hits on psegment):
988 
989  if( Psegments[ pseg_pos ].size() > 2 ) {
990 
991  // looks more exciting than it is. Here the weight is calculated. It is the difference in x of the last two and one but the last two hits,
992  // divided by the distance of the corresponding hits. Please refer to twiki page XXXX or CMS Note YYY (and use layer_distance)
993 
994  weight_A[ pseg_pos ] += theWeight(
995  (*(Psegments[ pseg_pos ].end()-1 ))->localPosition().x(),
996  (*(Psegments[ pseg_pos ].end()-2 ))->localPosition().x(),
997  (*(Psegments[ pseg_pos ].end()-3 ))->localPosition().x(),
998  float((*(Psegments[ pseg_pos ].end()-1))->cscDetId().layer()),
999  float((*(Psegments[ pseg_pos ].end()-2))->cscDetId().layer()),
1000  float((*(Psegments[ pseg_pos ].end()-3))->cscDetId().layer())
1001  );
1002 
1003  weight_B[ pseg_pos ] += theWeight(
1004  (*(Psegments[ pseg_pos ].end()-1 ))->localPosition().y(),
1005  (*(Psegments[ pseg_pos ].end()-2 ))->localPosition().y(),
1006  (*(Psegments[ pseg_pos ].end()-3 ))->localPosition().y(),
1007  float((*(Psegments[ pseg_pos ].end()-1))->cscDetId().layer()),
1008  float((*(Psegments[ pseg_pos ].end()-2))->cscDetId().layer()),
1009  float((*(Psegments[ pseg_pos ].end()-3))->cscDetId().layer())
1010  );
1011 
1012  // if we have picked up the last hit go looking for pseg with the lowest (and second lowest?) weight
1013 
1014  if(int(Psegments[ pseg_pos ].size()) == n_layers_occupied_tot) {
1015 
1016  curv_A[ pseg_pos ] += theWeight(
1017  (*(Psegments[ pseg_pos ].end()-1 ))->localPosition().x(),
1018  (*(Psegments[ pseg_pos ].end()-midlayer_pointer[n_layers_occupied_tot-1] ))->localPosition().x(),
1019  (*(Psegments[ pseg_pos ].end()-n_layers_occupied_tot ))->localPosition().x(),
1020  float((*(Psegments[ pseg_pos ].end()-1))->cscDetId().layer()),
1021  float((*(Psegments[ pseg_pos ].end()-midlayer_pointer[n_layers_occupied_tot-1] ))->cscDetId().layer()),
1022  float((*(Psegments[ pseg_pos ].end()-n_layers_occupied_tot ))->cscDetId().layer())
1023  );
1024 
1025  if (curv_A[ pseg_pos ] > curvePenaltyThreshold) weight_A[ pseg_pos ] = weight_A[ pseg_pos ] * curvePenalty;
1026 
1027  if (weight_B[ pseg_pos ] > a_yweightPenaltyThreshold[thestation][thering]) weight_A[ pseg_pos ] = weight_A[ pseg_pos ] * yweightPenalty;
1028 
1029  if (weight_A[ pseg_pos ] < min_weight_A ) {
1030  min_weight_A = weight_A[ pseg_pos ];
1031  //best_weight_B = weight_B[ pseg_pos ];
1032  //best_curv_A = curv_A[ pseg_pos ];
1033  best_pseg = pseg_pos ;
1034  }
1035 
1036  }
1037 
1038  // alternative: fill map with weight and pseg (which is already ordered)? Seems a very good tool to go looking for segments from.
1039  // As I understand, the segments would be inserted according to their weight, so the list would "automatically" be sorted.
1040 
1041  }
1042 
1043  if ( n_layers_occupied_tot > 3 ) {
1044  if (pseg < orig_number_of_noL1_psegs && (n_layers_processed != 2)) {
1045  if(( Psegments_noL1[ pseg_noL1_pos ].size() > 2 ) ) {
1046 
1047  // looks more exciting than it is. Here the weight is calculated. It is the difference in x of the last two and one but the last two hits,
1048  // divided by the distance of the corresponding hits. Please refer to twiki page XXXX or CMS Note YYY (and use layer_distance)
1049 
1050  weight_noL1_A[ pseg_noL1_pos ] += theWeight(
1051  (*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->localPosition().x(),
1052  (*(Psegments_noL1[ pseg_noL1_pos ].end()-2 ))->localPosition().x(),
1053  (*(Psegments_noL1[ pseg_noL1_pos ].end()-3 ))->localPosition().x(),
1054  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-1))->cscDetId().layer()),
1055  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-2))->cscDetId().layer()),
1056  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-3))->cscDetId().layer())
1057  );
1058 
1059  weight_noL1_B[ pseg_noL1_pos ] += theWeight(
1060  (*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->localPosition().y(),
1061  (*(Psegments_noL1[ pseg_noL1_pos ].end()-2 ))->localPosition().y(),
1062  (*(Psegments_noL1[ pseg_noL1_pos ].end()-3 ))->localPosition().y(),
1063  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-1))->cscDetId().layer()),
1064  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-2))->cscDetId().layer()),
1065  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-3))->cscDetId().layer())
1066  );
1067 
1068  //if we have picked up the last hit go looking for pseg with the lowest (and second lowest?) weight
1069 
1070  if(int(Psegments_noL1[ pseg_noL1_pos ].size()) == n_layers_occupied_tot -1 ) {
1071 
1072  curv_noL1_A[ pseg_noL1_pos ] += theWeight(
1073  (*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->localPosition().x(),
1074  (*(Psegments_noL1[ pseg_noL1_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
1075  (*(Psegments_noL1[ pseg_noL1_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
1076  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->cscDetId().layer()),
1077  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
1078  float((*(Psegments_noL1[ pseg_noL1_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
1079  );
1080 
1081  if (curv_noL1_A[ pseg_noL1_pos ] > curvePenaltyThreshold) weight_noL1_A[ pseg_noL1_pos ] = weight_noL1_A[ pseg_noL1_pos ] * curvePenalty;
1082 
1083  if (weight_noL1_B[ pseg_noL1_pos ] > a_yweightPenaltyThreshold[thestation][thering])
1084  weight_noL1_A[ pseg_noL1_pos ] = weight_noL1_A[ pseg_noL1_pos ] * yweightPenalty;
1085 
1086  if (weight_noL1_A[ pseg_noL1_pos ] < min_weight_noLx_A ) {
1087  min_weight_noLx_A = weight_noL1_A[ pseg_noL1_pos ];
1088  //best_weight_noLx_B = weight_noL1_B[ pseg_noL1_pos ];
1089  //best_curv_noLx_A = curv_noL1_A[ pseg_noL1_pos ];
1090  best_noLx_pseg = pseg_noL1_pos;
1091  best_Layer_noLx = 1;
1092  }
1093 
1094  }
1095 
1096  // alternative: fill map with weight and pseg (which is already ordered)? Seems a very good tool to go looking for segments from.
1097  // As I understand, the segments would be inserted according to their weight, so the list would "automatically" be sorted.
1098 
1099  }
1100  }
1101  }
1102 
1103  if ( n_layers_occupied_tot > 3 ) {
1104  if (pseg < orig_number_of_noL2_psegs && ( n_layers_processed != 2 )) {
1105  if(( Psegments_noL2[ pseg_noL2_pos ].size() > 2 )) {
1106 
1107  // looks more exciting than it is. Here the weight is calculated. It is the difference in x of the last two and one but the last two hits,
1108  // divided by the distance of the corresponding hits. Please refer to twiki page XXXX or CMS Note YYY (and use layer_distance)
1109 
1110  weight_noL2_A[ pseg_noL2_pos ] += theWeight(
1111  (*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->localPosition().x(),
1112  (*(Psegments_noL2[ pseg_noL2_pos ].end()-2 ))->localPosition().x(),
1113  (*(Psegments_noL2[ pseg_noL2_pos ].end()-3 ))->localPosition().x(),
1114  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-1))->cscDetId().layer()),
1115  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-2))->cscDetId().layer()),
1116  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-3))->cscDetId().layer())
1117  );
1118 
1119  weight_noL2_B[ pseg_noL2_pos ] += theWeight(
1120  (*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->localPosition().y(),
1121  (*(Psegments_noL2[ pseg_noL2_pos ].end()-2 ))->localPosition().y(),
1122  (*(Psegments_noL2[ pseg_noL2_pos ].end()-3 ))->localPosition().y(),
1123  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-1))->cscDetId().layer()),
1124  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-2))->cscDetId().layer()),
1125  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-3))->cscDetId().layer())
1126  );
1127 
1128  //if we have picked up the last hit go looking for pseg with the lowest (and second lowest?) weight
1129 
1130  if(int(Psegments_noL2[ pseg_noL2_pos ].size()) == n_layers_occupied_tot -1 ) {
1131 
1132  curv_noL2_A[ pseg_noL2_pos ] += theWeight(
1133  (*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->localPosition().x(),
1134  (*(Psegments_noL2[ pseg_noL2_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
1135  (*(Psegments_noL2[ pseg_noL2_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
1136  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->cscDetId().layer()),
1137  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
1138  float((*(Psegments_noL2[ pseg_noL2_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
1139  );
1140 
1141  if (curv_noL2_A[ pseg_noL2_pos ] > curvePenaltyThreshold) weight_noL2_A[ pseg_noL2_pos ] = weight_noL2_A[ pseg_noL2_pos ] * curvePenalty;
1142 
1143  if (weight_noL2_B[ pseg_noL2_pos ] > a_yweightPenaltyThreshold[thestation][thering])
1144  weight_noL2_A[ pseg_noL2_pos ] = weight_noL2_A[ pseg_noL2_pos ] * yweightPenalty;
1145 
1146  if (weight_noL2_A[ pseg_noL2_pos ] < min_weight_noLx_A ) {
1147  min_weight_noLx_A = weight_noL2_A[ pseg_noL2_pos ];
1148  //best_weight_noLx_B = weight_noL2_B[ pseg_noL2_pos ];
1149  //best_curv_noLx_A = curv_noL2_A[ pseg_noL2_pos ];
1150  best_noLx_pseg = pseg_noL2_pos;
1151  best_Layer_noLx = 2;
1152  }
1153 
1154  }
1155 
1156  // alternative: fill map with weight and pseg (which is already ordered)? Seems a very good tool to go looking for segments from.
1157  // As I understand, the segments would be inserted according to their weight, so the list would "automatically" be sorted.
1158 
1159  }
1160  }
1161  }
1162 
1163  if ( n_layers_occupied_tot > 3 ) {
1164  if (pseg < orig_number_of_noL3_psegs && ( n_layers_processed != 3 )) {
1165  if(( Psegments_noL3[ pseg_noL3_pos ].size() > 2 )) {
1166 
1167  // looks more exciting than it is. Here the weight is calculated. It is the difference in x of the last two and one but the last two hits,
1168  // divided by the distance of the corresponding hits. Please refer to twiki page XXXX or CMS Note YYY (and use layer_distance)
1169 
1170  weight_noL3_A[ pseg_noL3_pos ] += theWeight(
1171  (*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->localPosition().x(),
1172  (*(Psegments_noL3[ pseg_noL3_pos ].end()-2 ))->localPosition().x(),
1173  (*(Psegments_noL3[ pseg_noL3_pos ].end()-3 ))->localPosition().x(),
1174  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-1))->cscDetId().layer()),
1175  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-2))->cscDetId().layer()),
1176  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-3))->cscDetId().layer())
1177  );
1178 
1179  weight_noL3_B[ pseg_noL3_pos ] += theWeight(
1180  (*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->localPosition().y(),
1181  (*(Psegments_noL3[ pseg_noL3_pos ].end()-2 ))->localPosition().y(),
1182  (*(Psegments_noL3[ pseg_noL3_pos ].end()-3 ))->localPosition().y(),
1183  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-1))->cscDetId().layer()),
1184  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-2))->cscDetId().layer()),
1185  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-3))->cscDetId().layer())
1186  );
1187 
1188  //if we have picked up the last hit go looking for pseg with the lowest (and second lowest?) weight
1189 
1190  if(int(Psegments_noL3[ pseg_noL3_pos ].size()) == n_layers_occupied_tot -1 ) {
1191 
1192  curv_noL3_A[ pseg_noL3_pos ] += theWeight(
1193  (*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->localPosition().x(),
1194  (*(Psegments_noL3[ pseg_noL3_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
1195  (*(Psegments_noL3[ pseg_noL3_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
1196  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->cscDetId().layer()),
1197  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
1198  float((*(Psegments_noL3[ pseg_noL3_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
1199  );
1200 
1201  if (curv_noL3_A[ pseg_noL3_pos ] > curvePenaltyThreshold) weight_noL3_A[ pseg_noL3_pos ] = weight_noL3_A[ pseg_noL3_pos ] * curvePenalty;
1202 
1203  if (weight_noL3_B[ pseg_noL3_pos ] > a_yweightPenaltyThreshold[thestation][thering])
1204  weight_noL3_A[ pseg_noL3_pos ] = weight_noL3_A[ pseg_noL3_pos ] * yweightPenalty;
1205 
1206  if (weight_noL3_A[ pseg_noL3_pos ] < min_weight_noLx_A ) {
1207  min_weight_noLx_A = weight_noL3_A[ pseg_noL3_pos ];
1208  //best_weight_noLx_B = weight_noL3_B[ pseg_noL3_pos ];
1209  //best_curv_noLx_A = curv_noL3_A[ pseg_noL3_pos ];
1210  best_noLx_pseg = pseg_noL3_pos;
1211  best_Layer_noLx = 3;
1212  }
1213 
1214  }
1215 
1216  // alternative: fill map with weight and pseg (which is already ordered)? Seems a very good tool to go looking for segments from.
1217  // As I understand, the segments would be inserted according to their weight, so the list would "automatically" be sorted.
1218 
1219  }
1220  }
1221  }
1222 
1223  if ( n_layers_occupied_tot > 3 ) {
1224  if (pseg < orig_number_of_noL4_psegs && ( n_layers_processed != 4 )) {
1225  if(( Psegments_noL4[ pseg_noL4_pos ].size() > 2 )) {
1226 
1227  // looks more exciting than it is. Here the weight is calculated. It is the difference in x of the last two and one but the last two hits,
1228  // divided by the distance of the corresponding hits. Please refer to twiki page XXXX or CMS Note YYY (and use layer_distance)
1229 
1230  weight_noL4_A[ pseg_noL4_pos ] += theWeight(
1231  (*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->localPosition().x(),
1232  (*(Psegments_noL4[ pseg_noL4_pos ].end()-2 ))->localPosition().x(),
1233  (*(Psegments_noL4[ pseg_noL4_pos ].end()-3 ))->localPosition().x(),
1234  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-1))->cscDetId().layer()),
1235  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-2))->cscDetId().layer()),
1236  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-3))->cscDetId().layer())
1237  );
1238 
1239  weight_noL4_B[ pseg_noL4_pos ] += theWeight(
1240  (*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->localPosition().y(),
1241  (*(Psegments_noL4[ pseg_noL4_pos ].end()-2 ))->localPosition().y(),
1242  (*(Psegments_noL4[ pseg_noL4_pos ].end()-3 ))->localPosition().y(),
1243  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-1))->cscDetId().layer()),
1244  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-2))->cscDetId().layer()),
1245  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-3))->cscDetId().layer())
1246  );
1247 
1248  //if we have picked up the last hit go looking for pseg with the lowest (and second lowest?) weight
1249 
1250  if(int(Psegments_noL4[ pseg_noL4_pos ].size()) == n_layers_occupied_tot -1 ) {
1251 
1252  curv_noL4_A[ pseg_noL4_pos ] += theWeight(
1253  (*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->localPosition().x(),
1254  (*(Psegments_noL4[ pseg_noL4_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
1255  (*(Psegments_noL4[ pseg_noL4_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
1256  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->cscDetId().layer()),
1257  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
1258  float((*(Psegments_noL4[ pseg_noL4_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
1259  );
1260 
1261  if (curv_noL4_A[ pseg_noL4_pos ] > curvePenaltyThreshold) weight_noL4_A[ pseg_noL4_pos ] = weight_noL4_A[ pseg_noL4_pos ] * curvePenalty;
1262 
1263  if (weight_noL4_B[ pseg_noL4_pos ] > a_yweightPenaltyThreshold[thestation][thering])
1264  weight_noL4_A[ pseg_noL4_pos ] = weight_noL4_A[ pseg_noL4_pos ] * yweightPenalty;
1265 
1266  if (weight_noL4_A[ pseg_noL4_pos ] < min_weight_noLx_A ) {
1267  min_weight_noLx_A = weight_noL4_A[ pseg_noL4_pos ];
1268  //best_weight_noLx_B = weight_noL4_B[ pseg_noL4_pos ];
1269  //best_curv_noLx_A = curv_noL4_A[ pseg_noL4_pos ];
1270  best_noLx_pseg = pseg_noL4_pos;
1271  best_Layer_noLx = 4;
1272  }
1273 
1274  }
1275 
1276  // alternative: fill map with weight and pseg (which is already ordered)? Seems a very good tool to go looking for segments from.
1277  // As I understand, the segments would be inserted according to their weight, so the list would "automatically" be sorted.
1278 
1279  }
1280  }
1281  }
1282 
1283  if ( n_layers_occupied_tot > 4 ) {
1284  if (pseg < orig_number_of_noL5_psegs && ( n_layers_processed != 5 )) {
1285  if(( Psegments_noL5[ pseg_noL5_pos ].size() > 2 )){
1286 
1287  // looks more exciting than it is. Here the weight is calculated. It is the difference in x of the last two and one but the last two hits,
1288  // divided by the distance of the corresponding hits. Please refer to twiki page XXXX or CMS Note YYY (and use layer_distance)
1289 
1290  weight_noL5_A[ pseg_noL5_pos ] += theWeight(
1291  (*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->localPosition().x(),
1292  (*(Psegments_noL5[ pseg_noL5_pos ].end()-2 ))->localPosition().x(),
1293  (*(Psegments_noL5[ pseg_noL5_pos ].end()-3 ))->localPosition().x(),
1294  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-1))->cscDetId().layer()),
1295  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-2))->cscDetId().layer()),
1296  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-3))->cscDetId().layer())
1297  );
1298 
1299  weight_noL5_B[ pseg_noL5_pos ] += theWeight(
1300  (*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->localPosition().y(),
1301  (*(Psegments_noL5[ pseg_noL5_pos ].end()-2 ))->localPosition().y(),
1302  (*(Psegments_noL5[ pseg_noL5_pos ].end()-3 ))->localPosition().y(),
1303  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-1))->cscDetId().layer()),
1304  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-2))->cscDetId().layer()),
1305  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-3))->cscDetId().layer())
1306  );
1307 
1308  //if we have picked up the last hit go looking for pseg with the lowest (and second lowest?) weight
1309 
1310  if(int(Psegments_noL5[ pseg_noL5_pos ].size()) == n_layers_occupied_tot -1 ) {
1311 
1312  curv_noL5_A[ pseg_noL5_pos ] += theWeight(
1313  (*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->localPosition().x(),
1314  (*(Psegments_noL5[ pseg_noL5_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
1315  (*(Psegments_noL5[ pseg_noL5_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
1316  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->cscDetId().layer()),
1317  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
1318  float((*(Psegments_noL5[ pseg_noL5_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
1319  );
1320 
1321  if (curv_noL5_A[ pseg_noL5_pos ] > curvePenaltyThreshold) weight_noL5_A[ pseg_noL5_pos ] = weight_noL5_A[ pseg_noL5_pos ] * curvePenalty;
1322 
1323  if (weight_noL5_B[ pseg_noL5_pos ] > a_yweightPenaltyThreshold[thestation][thering])
1324  weight_noL5_A[ pseg_noL5_pos ] = weight_noL5_A[ pseg_noL5_pos ] * yweightPenalty;
1325 
1326  if (weight_noL5_A[ pseg_noL5_pos ] < min_weight_noLx_A ) {
1327  min_weight_noLx_A = weight_noL5_A[ pseg_noL5_pos ];
1328  //best_weight_noLx_B = weight_noL5_B[ pseg_noL5_pos ];
1329  //best_curv_noLx_A = curv_noL5_A[ pseg_noL5_pos ];
1330  best_noLx_pseg = pseg_noL5_pos;
1331  best_Layer_noLx = 5;
1332  }
1333 
1334  }
1335 
1336  // alternative: fill map with weight and pseg (which is already ordered)? Seems a very good tool to go looking for segments from.
1337  // As I understand, the segments would be inserted according to their weight, so the list would "automatically" be sorted.
1338 
1339  }
1340  }
1341  }
1342 
1343  if ( n_layers_occupied_tot > 5 ) {
1344  if (pseg < orig_number_of_noL6_psegs && ( n_layers_processed != 6 )) {
1345  if(( Psegments_noL6[ pseg_noL6_pos ].size() > 2 )){
1346 
1347  // looks more exciting than it is. Here the weight is calculated. It is the difference in x of the last two and one but the last two hits,
1348  // divided by the distance of the corresponding hits. Please refer to twiki page XXXX or CMS Note YYY (and use layer_distance)
1349 
1350  weight_noL6_A[ pseg_noL6_pos ] += theWeight(
1351  (*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->localPosition().x(),
1352  (*(Psegments_noL6[ pseg_noL6_pos ].end()-2 ))->localPosition().x(),
1353  (*(Psegments_noL6[ pseg_noL6_pos ].end()-3 ))->localPosition().x(),
1354  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-1))->cscDetId().layer()),
1355  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-2))->cscDetId().layer()),
1356  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-3))->cscDetId().layer())
1357  );
1358 
1359  weight_noL6_B[ pseg_noL6_pos ] += theWeight(
1360  (*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->localPosition().y(),
1361  (*(Psegments_noL6[ pseg_noL6_pos ].end()-2 ))->localPosition().y(),
1362  (*(Psegments_noL6[ pseg_noL6_pos ].end()-3 ))->localPosition().y(),
1363  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-1))->cscDetId().layer()),
1364  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-2))->cscDetId().layer()),
1365  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-3))->cscDetId().layer())
1366  );
1367 
1368  //if we have picked up the last hit go looking for pseg with the lowest (and second lowest?) weight
1369 
1370  if(int(Psegments_noL6[ pseg_noL6_pos ].size()) == n_layers_occupied_tot -1 ) {
1371 
1372  curv_noL6_A[ pseg_noL6_pos ] += theWeight(
1373  (*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->localPosition().x(),
1374  (*(Psegments_noL6[ pseg_noL6_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
1375  (*(Psegments_noL6[ pseg_noL6_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
1376  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->cscDetId().layer()),
1377  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
1378  float((*(Psegments_noL6[ pseg_noL6_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
1379  );
1380 
1381  if (curv_noL6_A[ pseg_noL6_pos ] > curvePenaltyThreshold) weight_noL6_A[ pseg_noL6_pos ] = weight_noL6_A[ pseg_noL6_pos ] * curvePenalty;
1382 
1383  if (weight_noL6_B[ pseg_noL6_pos ] > a_yweightPenaltyThreshold[thestation][thering])
1384  weight_noL6_A[ pseg_noL6_pos ] = weight_noL6_A[ pseg_noL6_pos ] * yweightPenalty;
1385 
1386  if (weight_noL6_A[ pseg_noL6_pos ] < min_weight_noLx_A ) {
1387  min_weight_noLx_A = weight_noL6_A[ pseg_noL6_pos ];
1388  //best_weight_noLx_B = weight_noL6_B[ pseg_noL6_pos ];
1389  //best_curv_noLx_A = curv_noL6_A[ pseg_noL6_pos ];
1390  best_noLx_pseg = pseg_noL6_pos;
1391  best_Layer_noLx = 6;
1392  }
1393 
1394  }
1395 
1396  // alternative: fill map with weight and pseg (which is already ordered)? Seems a very good tool to go looking for segments from.
1397  // As I understand, the segments would be inserted according to their weight, so the list would "automatically" be sorted.
1398 
1399  }
1400  }
1401  }
1402 
1403  }
1404  }
1405  }
1406  }
1407 
1408  //************************************************************************;
1409  //*** End segment building *******************************************;
1410  //************************************************************************;
1411 
1412  // Important part! Here segment(s) are actually chosen. All the good segments
1413  // could be chosen or some (best) ones only (in order to save time).
1414 
1415  // Check if there is a segment with n-1 hits that has a signifcantly better
1416  // weight than the best n hit segment
1417 
1418  // IBL 070828: implicit assumption here is that there will always only be one segment per
1419  // cluster - if there are >1 we will need to find out which segment the alternative n-1 hit
1420  // protosegment belongs to!
1421 
1422 
1423  //float chosen_weight = min_weight_A;
1424  //float chosen_ywgt = best_weight_B;
1425  //float chosen_curv = best_curv_A;
1426  //int chosen_nlayers = n_layers_occupied_tot;
1427  int chosen_pseg = best_pseg;
1428  if (best_pseg<0) {
1429  return segmentInChamber;
1430  }
1433 
1434  float hit_drop_limit = -999999.999;
1435 
1436  // define different weight improvement requirements depending on how many layers are in the segment candidate
1437  switch ( n_layers_processed ) {
1438  case 1 :
1439  // do nothing;
1440  break;
1441  case 2 :
1442  // do nothing;
1443  break;
1444  case 3 :
1445  // do nothing;
1446  break;
1447  case 4 :
1448  hit_drop_limit = hitDropLimit6Hits * (1./2.) * hitDropLimit4Hits;
1449  if ((best_Layer_noLx < 1) || (best_Layer_noLx > 4)) {
1450  // std::cout<<"CSCSegAlgoST: For four layers, best_Layer_noLx = "<< best_Layer_noLx << std::endl;
1451  }
1452  if ((best_Layer_noLx == 2) || (best_Layer_noLx == 3)) hit_drop_limit = hit_drop_limit * (1./2.);
1453  break;
1454  case 5 :
1455  hit_drop_limit = hitDropLimit6Hits * (2./3.) * hitDropLimit5Hits;
1456  if ((best_Layer_noLx < 1) || (best_Layer_noLx > 5)) {
1457  // std::cout<<"CSCSegAlgoST: For five layers, best_Layer_noLx = "<< best_Layer_noLx << std::endl;
1458  }
1459  if ((best_Layer_noLx == 2) || (best_Layer_noLx == 4)) hit_drop_limit = hit_drop_limit * (1./2.);
1460  if (best_Layer_noLx == 3) hit_drop_limit = hit_drop_limit * (1./3.);
1461  break;
1462  case 6 :
1463  hit_drop_limit = hitDropLimit6Hits * (3./4.);
1464  if ((best_Layer_noLx < 1) || (best_Layer_noLx > 6)) {
1465  // std::cout<<"CSCSegAlgoST: For six layers, best_Layer_noLx = "<< best_Layer_noLx << std::endl;
1466  }
1467  if ((best_Layer_noLx == 2) || (best_Layer_noLx == 5)) hit_drop_limit = hit_drop_limit * (1./2.);
1468  if ((best_Layer_noLx == 3) || (best_Layer_noLx == 4)) hit_drop_limit = hit_drop_limit * (1./3.);
1469  break;
1470 
1471  default :
1472  // Fallback - should never occur.
1473  LogDebug("CSCSegment|CSC") <<"CSCSegAlgoST: Unexpected number of layers with hits - please inform developers.";
1474  // std::cout<<"CSCSegAlgoST: Unexpected number of layers with hits - please inform developers."<<std::endl;
1475  hit_drop_limit = 0.1;
1476  }
1477 
1478  // choose the NoLx collection (the one that contains the best N-1 candidate)
1479  switch ( best_Layer_noLx ) {
1480  case 1 :
1481  Psegments_noLx.clear();
1483  weight_noLx_A.clear();
1485  break;
1486  case 2 :
1487  Psegments_noLx.clear();
1489  weight_noLx_A.clear();
1491  break;
1492  case 3 :
1493  Psegments_noLx.clear();
1495  weight_noLx_A.clear();
1497  break;
1498  case 4 :
1499  Psegments_noLx.clear();
1501  weight_noLx_A.clear();
1503  break;
1504  case 5 :
1505  Psegments_noLx.clear();
1507  weight_noLx_A.clear();
1509  break;
1510  case 6 :
1511  Psegments_noLx.clear();
1513  weight_noLx_A.clear();
1515  break;
1516 
1517  default :
1518  // Fallback - should occur only for preclusters with only 3 layers with hits.
1519  Psegments_noLx.clear();
1520  weight_noLx_A.clear();
1521  }
1522 
1523  if( min_weight_A > 0. ) {
1524  if ( min_weight_noLx_A/min_weight_A < hit_drop_limit ) {
1525  //chosen_weight = min_weight_noLx_A;
1526  //chosen_ywgt = best_weight_noLx_B;
1527  //chosen_curv = best_curv_noLx_A;
1528  //chosen_nlayers = n_layers_occupied_tot-1;
1529  chosen_pseg = best_noLx_pseg;
1530  chosen_Psegments.clear();
1531  chosen_weight_A.clear();
1534  }
1535  }
1536 
1537  if(onlyBestSegment) {
1538  ChooseSegments2a( chosen_Psegments, chosen_pseg );
1539  }
1540  else {
1542  }
1543 
1544  for(unsigned int iSegment=0; iSegment<GoodSegments.size();++iSegment){
1545  protoSegment = GoodSegments[iSegment];
1546  passCondNumber=false;
1547  passCondNumber_2 = false;
1548  protoChiUCorrection=1.0;
1549  doSlopesAndChi2();
1550  // Attempt to handle numerical instability of the fit;
1551  // Any segment with protoChi2/protoNDF>chi2Norm_3D_
1552  // considered as that one potentially suffering from
1553  // numerical instability in fit.
1554  if(correctCov_){
1555  // Call the fit with prefitting option;
1556  // First fit a straight line to X-Z coordinates
1557  // and calculate chi^2 (chiUZ in correctTheCovX(void)) for X-Z fit;
1558  // Scale up errors in X if chiUZ too big (default 20);
1559  // Refit XY-Z with the scaled up X errors
1561  passCondNumber = true;
1562  doSlopesAndChi2();
1563  }
1564  if(protoChiUCorrection<1.00005){
1565  LogDebug("CSCSegment|segmWierd") << "Wierd segment, ErrXX scaled, refit " <<std::endl;
1567  // Call the fit with direct adjustment of condition number;
1568  // If the attempt to improve fit by scaling up X error fails
1569  // call the procedure to make the condition number of M compatible with
1570  // the precision of X and Y measurements;
1571  // Achieved by decreasing abs value of the Covariance
1572  LogDebug("CSCSegment|segmWierd") << "Wierd segment, ErrXY changed to match cond. number, refit " << std::endl;
1573  passCondNumber_2=true;
1574  doSlopesAndChi2();
1575  }
1576  }
1577  // Call the pre-pruning procedure;
1578  // If the attempt to improve fit by scaling up X error is successfull,
1579  // while scale factor for X errors is too big.
1580  // Prune the recHit inducing the biggest contribution into X-Z chi^2
1581  // and refit;
1583  (protoSegment.size()>3)){
1584  LogDebug("CSCSegment|segmWierd") << "Scale factor protoChiUCorrection too big, pre-Prune, refit " << std::endl;
1585  protoSegment.erase(protoSegment.begin()+(maxContrIndex),
1586  protoSegment.begin()+(maxContrIndex+1));
1587  doSlopesAndChi2();
1588  }
1589  }
1590 
1592  // calculate error matrix
1593  AlgebraicSymMatrix protoErrors = calculateError();
1594  // but reorder components to match what's required by TrackingRecHit interface
1595  // i.e. slopes first, then positions
1596  flipErrors( protoErrors );
1597  //
1599  segmentInChamber.push_back(temp);
1600  }
1601  return segmentInChamber;
1602 }
1603 
1604 void CSCSegAlgoST::ChooseSegments2a(std::vector< ChamberHitContainer > & chosen_segments, int chosen_seg) {
1605  // just return best segment
1606  GoodSegments.clear();
1607  GoodSegments.push_back( chosen_segments[chosen_seg] );
1608 }
1609 
1610 void CSCSegAlgoST::ChooseSegments3(std::vector< ChamberHitContainer > & chosen_segments, std::vector< float > & chosen_weight, int chosen_seg) {
1611 
1612  int SumCommonHits = 0;
1613  GoodSegments.clear();
1614  int nr_remaining_candidates;
1615  unsigned int nr_of_segment_candidates;
1616 
1617  nr_remaining_candidates = nr_of_segment_candidates = chosen_segments.size();
1618 
1619  // always select and return best protosegment:
1620  GoodSegments.push_back( chosen_segments[ chosen_seg ] );
1621 
1622  float chosen_weight_temp = 999999.;
1623  int chosen_seg_temp = -1;
1624 
1625  // try to find further segment candidates:
1626  while( nr_remaining_candidates > 0 ) {
1627 
1628  for(unsigned int iCand=0; iCand < nr_of_segment_candidates; ++iCand) {
1629  //only compare current best to psegs that have not been marked bad:
1630  if( chosen_weight[iCand] < 0. ) continue;
1631  SumCommonHits = 0;
1632 
1633  for( int ihits = 0; ihits < int(chosen_segments[iCand].size()); ++ihits ) { // iCand and iiCand NEED to have same nr of hits! (always have by construction)
1634  if( chosen_segments[iCand][ihits] == chosen_segments[chosen_seg][ihits]) {
1635  ++SumCommonHits;
1636  }
1637  }
1638 
1639  //mark a pseg bad:
1640  if(SumCommonHits>1) { // needs to be a card; should be investigated first
1641  chosen_weight[iCand] = -1.;
1642  nr_remaining_candidates -= 1;
1643  }
1644  else {
1645  // save the protosegment with the smallest weight
1646  if( chosen_weight[ iCand ] < chosen_weight_temp ) {
1647  chosen_weight_temp = chosen_weight[ iCand ];
1648  chosen_seg_temp = iCand ;
1649  }
1650  }
1651  }
1652 
1653  if( chosen_seg_temp > -1 ) GoodSegments.push_back( chosen_segments[ chosen_seg_temp ] );
1654 
1655  chosen_seg = chosen_seg_temp;
1656  // re-initialze temporary best parameters
1657  chosen_weight_temp = 999999;
1658  chosen_seg_temp = -1;
1659  }
1660 }
1661 
1662 void CSCSegAlgoST::ChooseSegments2(int best_seg) {
1663  // std::vector <int> CommonHits(6); // nice concept :)
1664  std::vector <unsigned int> BadCandidate;
1665  int SumCommonHits =0;
1666  GoodSegments.clear();
1667  BadCandidate.clear();
1668  for(unsigned int iCand=0;iCand<Psegments.size();++iCand) {
1669  // skip here if segment was marked bad
1670  for(unsigned int iiCand=iCand+1;iiCand<Psegments.size();++iiCand){
1671  // skip here too if segment was marked bad
1672  SumCommonHits =0;
1673  if( Psegments[iCand].size() != Psegments[iiCand].size() ) {
1674  LogDebug("CSCSegment|CSC") <<"CSCSegmentST::ChooseSegments2: ALARM!! THIS should not happen!!";
1675 // std::cout<<"CSCSegmentST::ChooseSegments2: ALARM!! THIS should not happen!!"<<std::endl;
1676  }
1677  else {
1678  for( int ihits = 0; ihits < int(Psegments[iCand].size()); ++ihits ) { // iCand and iiCand NEED to have same nr of hits! (alsways have by construction)
1679  if( Psegments[iCand][ihits] == Psegments[iiCand][ihits]) {
1680  ++SumCommonHits;
1681  }
1682  }
1683  }
1684  if(SumCommonHits>1) {
1685  if( weight_A[iCand]>weight_A[iiCand] ) { // use weight_A here
1686  BadCandidate.push_back(iCand);
1687  // rather mark segment bad by an array which is in sync with protosegments!! e.g. set weight = weight*1000 or have an addidional array or set it to weight *= -1
1688  }
1689  else{
1690  BadCandidate.push_back(iiCand);
1691  // rather mark segment bad by an array which is in sync with protosegments!! e.g. set weight = weight*1000 or have an addidional array or set it to weight *= -1
1692  }
1693  }
1694  }
1695  }
1696  bool discard;
1697  for(unsigned int isegm=0;isegm<Psegments.size();++isegm) {
1698  // For best results another iteration/comparison over Psegments
1699  //should be applied here... It would make the program much slower.
1700  discard = false;
1701  for(unsigned int ibad=0;ibad<BadCandidate.size();++ibad) {
1702  // can save this loop if we used an array in sync with Psegments!!!!
1703  if(isegm == BadCandidate[ibad]) {
1704  discard = true;
1705  }
1706  }
1707  if(!discard) {
1708  GoodSegments.push_back( Psegments[isegm] );
1709  }
1710  }
1711 }
1712 //Method doSlopesAndChi2
1713 // fitSlopes() and fillChiSquared() are always called one after the other
1714 // In fact the code is duplicated in the two functions (as we need 2 loops) -
1715 // it is much better to fix that at some point
1717  fitSlopes();
1718  fillChiSquared();
1719 }
1720 /* Method fitSlopes
1721  *
1722  * Perform a Least Square Fit on a segment as per SK algo
1723  *
1724  */
1726  e_Cxx.clear();
1728  correctTheCovX();
1729  if(e_Cxx.size()!=protoSegment.size()){
1730  LogDebug("CSCSegment|segmWierd") << "e_Cxx.size()!=protoSegment.size() IT IS A SERIOUS PROBLEM!!! " <<std::endl;
1731  }
1732  }
1733  CLHEP::HepMatrix M(4,4,0);
1734  CLHEP::HepVector B(4,0);
1735  ChamberHitContainer::const_iterator ih = protoSegment.begin();
1736  for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
1737  const CSCRecHit2D& hit = (**ih);
1738  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
1739  GlobalPoint gp = layer->toGlobal(hit.localPosition());
1740  LocalPoint lp = theChamber->toLocal(gp);
1741  // ptc: Local position of hit w.r.t. chamber
1742  double u = lp.x();
1743  double v = lp.y();
1744  double z = lp.z();
1745  // ptc: Covariance matrix of local errors
1746  CLHEP::HepMatrix IC(2,2);
1748  IC(1,1) = e_Cxx.at(ih-protoSegment.begin());
1749  }
1750  else{
1751  IC(1,1) = hit.localPositionError().xx();
1752  }
1753  // IC(1,1) = hit.localPositionError().xx();
1754  IC(1,2) = hit.localPositionError().xy();
1755  IC(2,2) = hit.localPositionError().yy();
1756  IC(2,1) = IC(1,2); // since Cov is symmetric
1758  if(passCondNumber_2){
1759  correctTheCovMatrix(IC);
1760  }
1761  // ptc: Invert covariance matrix (and trap if it fails!)
1762  int ierr = 0;
1763  IC.invert(ierr); // inverts in place
1764  if (ierr != 0) {
1765  LogDebug("CSCSegment|CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC;
1766  // std::cout<< "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n"<<std::endl;
1767  }
1768 
1769  M(1,1) += IC(1,1);
1770  M(1,2) += IC(1,2);
1771  M(1,3) += IC(1,1) * z;
1772  M(1,4) += IC(1,2) * z;
1773  B(1) += u * IC(1,1) + v * IC(1,2);
1774 
1775  M(2,1) += IC(2,1);
1776  M(2,2) += IC(2,2);
1777  M(2,3) += IC(2,1) * z;
1778  M(2,4) += IC(2,2) * z;
1779  B(2) += u * IC(2,1) + v * IC(2,2);
1780 
1781  M(3,1) += IC(1,1) * z;
1782  M(3,2) += IC(1,2) * z;
1783  M(3,3) += IC(1,1) * z * z;
1784  M(3,4) += IC(1,2) * z * z;
1785  B(3) += ( u * IC(1,1) + v * IC(1,2) ) * z;
1786 
1787  M(4,1) += IC(2,1) * z;
1788  M(4,2) += IC(2,2) * z;
1789  M(4,3) += IC(2,1) * z * z;
1790  M(4,4) += IC(2,2) * z * z;
1791  B(4) += ( u * IC(2,1) + v * IC(2,2) ) * z;
1792  }
1793  CLHEP::HepVector p = solve(M, B);
1794 
1795  // Update member variables
1796  // Note that origin has local z = 0
1797  protoIntercept = LocalPoint(p(1), p(2), 0.);
1798  protoSlope_u = p(3);
1799  protoSlope_v = p(4);
1800 }
1801 /* Method fillChiSquared
1802  *
1803  * Determine Chi^2 for the proto wire segment
1804  *
1805  */
1807 
1808  double chsq = 0.;
1809 
1810  ChamberHitContainer::const_iterator ih;
1811  for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
1812 
1813  const CSCRecHit2D& hit = (**ih);
1814  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
1815  GlobalPoint gp = layer->toGlobal(hit.localPosition());
1816  LocalPoint lp = theChamber->toLocal(gp);
1817 
1818  double u = lp.x();
1819  double v = lp.y();
1820  double z = lp.z();
1821 
1822  double du = protoIntercept.x() + protoSlope_u * z - u;
1823  double dv = protoIntercept.y() + protoSlope_v * z - v;
1824 
1825  CLHEP::HepMatrix IC(2,2);
1827  IC(1,1) = e_Cxx.at(ih-protoSegment.begin());
1828  }
1829  else{
1830  IC(1,1) = hit.localPositionError().xx();
1831  }
1832  // IC(1,1) = hit.localPositionError().xx();
1833  IC(1,2) = hit.localPositionError().xy();
1834  IC(2,2) = hit.localPositionError().yy();
1835  IC(2,1) = IC(1,2);
1837  if(passCondNumber_2){
1838  correctTheCovMatrix(IC);
1839  }
1840 
1841  // Invert covariance matrix
1842  int ierr = 0;
1843  IC.invert(ierr);
1844  if (ierr != 0) {
1845  LogDebug("CSCSegment|CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC;
1846  // std::cout << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";
1847 
1848  }
1849 
1850  chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
1851  }
1852 
1853  protoChi2 = chsq;
1854  protoNDF = 2.*protoSegment.size() - 4;
1855 }
1856 /* fillLocalDirection
1857  *
1858  */
1860  // Always enforce direction of segment to point from IP outwards
1861  // (Incorrect for particles not coming from IP, of course.)
1862 
1863  double dxdz = protoSlope_u;
1864  double dydz = protoSlope_v;
1865  double dz = 1./sqrt(1. + dxdz*dxdz + dydz*dydz);
1866  double dx = dz*dxdz;
1867  double dy = dz*dydz;
1868  LocalVector localDir(dx,dy,dz);
1869 
1870  // localDir may need sign flip to ensure it points outward from IP
1871  // ptc: Examine its direction and origin in global z: to point outward
1872  // the localDir should always have same sign as global z...
1873 
1874  double globalZpos = ( theChamber->toGlobal( protoIntercept ) ).z();
1875  double globalZdir = ( theChamber->toGlobal( localDir ) ).z();
1876  double directionSign = globalZpos * globalZdir;
1877  protoDirection = (directionSign * localDir).unit();
1878 }
1879 /* weightMatrix
1880  *
1881  */
1883 
1884  std::vector<const CSCRecHit2D*>::const_iterator it;
1885  int nhits = protoSegment.size();
1886  AlgebraicSymMatrix matrix(2*nhits, 0);
1887  int row = 0;
1888 
1889  for (it = protoSegment.begin(); it != protoSegment.end(); ++it) {
1890 
1891  const CSCRecHit2D& hit = (**it);
1892  ++row;
1893  matrix(row, row) = protoChiUCorrection*hit.localPositionError().xx();
1894  matrix(row, row+1) = hit.localPositionError().xy();
1895  ++row;
1896  matrix(row, row-1) = hit.localPositionError().xy();
1897  matrix(row, row) = hit.localPositionError().yy();
1898  }
1899  int ierr;
1900  matrix.invert(ierr);
1901  return matrix;
1902 }
1903 
1904 
1905 /* derivativeMatrix
1906  *
1907  */
1908 CLHEP::HepMatrix CSCSegAlgoST::derivativeMatrix() const {
1909 
1910  ChamberHitContainer::const_iterator it;
1911  int nhits = protoSegment.size();
1912  CLHEP::HepMatrix matrix(2*nhits, 4);
1913  int row = 0;
1914 
1915  for(it = protoSegment.begin(); it != protoSegment.end(); ++it) {
1916 
1917  const CSCRecHit2D& hit = (**it);
1918  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
1919  GlobalPoint gp = layer->toGlobal(hit.localPosition());
1920  LocalPoint lp = theChamber->toLocal(gp);
1921  float z = lp.z();
1922  ++row;
1923  matrix(row, 1) = 1.;
1924  matrix(row, 3) = z;
1925  ++row;
1926  matrix(row, 2) = 1.;
1927  matrix(row, 4) = z;
1928  }
1929  return matrix;
1930 }
1931 
1932 
1933 /* calculateError
1934  *
1935  */
1937 
1940 
1941  // (AT W A)^-1
1942  // from http://www.phys.ufl.edu/~avery/fitting.html, part I
1943  int ierr;
1944  AlgebraicSymMatrix result = weights.similarityT(A);
1945  result.invert(ierr);
1946 
1947  // blithely assuming the inverting never fails...
1948  return result;
1949 }
1950 
1951 
1953 
1954  // The CSCSegment needs the error matrix re-arranged to match
1955  // parameters in order (uz, vz, u0, v0) where uz, vz = slopes, u0, v0 = intercepts
1956 
1957  AlgebraicSymMatrix hold( a );
1958 
1959  // errors on slopes into upper left
1960  a(1,1) = hold(3,3);
1961  a(1,2) = hold(3,4);
1962  a(2,1) = hold(4,3);
1963  a(2,2) = hold(4,4);
1964 
1965  // errors on positions into lower right
1966  a(3,3) = hold(1,1);
1967  a(3,4) = hold(1,2);
1968  a(4,3) = hold(2,1);
1969  a(4,4) = hold(2,2);
1970 
1971  // must also interchange off-diagonal elements of off-diagonal 2x2 submatrices
1972  a(4,1) = hold(2,3);
1973  a(3,2) = hold(1,4);
1974  a(2,3) = hold(4,1); // = hold(1,4)
1975  a(1,4) = hold(3,2); // = hold(2,3)
1976 }
1977 //
1979  std::vector<double> uu, vv, zz;
1980  //std::vector<double> e_Cxx;
1981  e_Cxx.clear();
1982  double sum_U_err=0.0;
1983  double sum_Z_U_err=0.0;
1984  double sum_Z2_U_err=0.0;
1985  double sum_U_U_err=0.0;
1986  double sum_UZ_U_err=0.0;
1987  std::vector<double> chiUZind;
1988  std::vector<double>::iterator chiContribution;
1989  double chiUZ=0.0;
1990  ChamberHitContainer::const_iterator ih = protoSegment.begin();
1991  for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
1992  const CSCRecHit2D& hit = (**ih);
1993  e_Cxx.push_back(hit.localPositionError().xx());
1994  //
1995  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
1996  GlobalPoint gp = layer->toGlobal(hit.localPosition());
1997  LocalPoint lp = theChamber->toLocal(gp);
1998  // ptc: Local position of hit w.r.t. chamber
1999  double u = lp.x();
2000  double v = lp.y();
2001  double z = lp.z();
2002  uu.push_back(u);
2003  vv.push_back(v);
2004  zz.push_back(z);
2006  sum_U_err += 1./e_Cxx.back();
2007  sum_Z_U_err += z/e_Cxx.back();
2008  sum_Z2_U_err += (z*z)/e_Cxx.back();
2009  sum_U_U_err += u/e_Cxx.back();
2010  sum_UZ_U_err += (u*z)/e_Cxx.back();
2011  }
2012 
2015 
2016  double denom=sum_U_err*sum_Z2_U_err-pow(sum_Z_U_err,2);
2017  double U0=(sum_Z2_U_err*sum_U_U_err-sum_Z_U_err*sum_UZ_U_err)/denom;
2018  double UZ=(sum_U_err*sum_UZ_U_err-sum_Z_U_err*sum_U_U_err)/denom;
2019 
2022 
2023  for(unsigned i=0; i<uu.size(); ++i){
2024  double uMean = U0+UZ*zz[i];
2025  chiUZind.push_back((pow((uMean-uu[i]),2))/e_Cxx[i]);
2026  chiUZ += (pow((uMean-uu[i]),2))/e_Cxx[i];
2027  }
2028  chiUZ = chiUZ/(uu.size()-2);
2029 
2030  if(chiUZ>=chi2Norm_2D_){
2032  for(unsigned i=0; i<uu.size(); ++i)
2034  }
2035 
2037 
2039  chiContribution=max_element(chiUZind.begin(),chiUZind.end());
2040  maxContrIndex = chiContribution - chiUZind.begin();
2041  /*
2042  for(unsigned i=0; i<chiUZind.size();++i){
2043  if(*chiContribution==chiUZind[i]){
2044  maxContrIndex=i;
2045  }
2046  }
2047  */
2048  }
2049  //
2050  //return e_Cxx;
2051 }
2052 //
2053 void CSCSegAlgoST::correctTheCovMatrix(CLHEP::HepMatrix &IC){
2054  //double condNumberCorr1=0.0;
2055  double condNumberCorr2=0.0;
2056  double detCov=0.0;
2057  double diag1=0.0;
2058  double diag2=0.0;
2059  double IC_12_corr=0.0;
2060  double IC_11_corr=0.0;
2061  if(!covToAnyNumberAll_){
2062  //condNumberCorr1=condSeed1_*IC(2,2);
2063  condNumberCorr2=condSeed2_*IC(2,2);
2064  diag1=IC(1,1)*IC(2,2);
2065  diag2=IC(1,2)*IC(1,2);
2066  detCov=fabs(diag1-diag2);
2067  if((diag1<condNumberCorr2)&&(diag2<condNumberCorr2)){
2068  if(covToAnyNumber_)
2069  IC(1,2)=covAnyNumber_;
2070  else{
2071  IC_11_corr=condSeed1_+fabs(IC(1,2))/IC(2,2);
2072  IC(1,1)=IC_11_corr;
2073  }
2074  }
2075 
2076  if(((detCov<condNumberCorr2)&&(diag1>condNumberCorr2))||
2077  ((diag2>condNumberCorr2)&&(detCov<condNumberCorr2)
2078  )){
2079  if(covToAnyNumber_)
2080  IC(1,2)=covAnyNumber_;
2081  else{
2082  IC_12_corr=sqrt(fabs(diag1-condNumberCorr2));
2083  if(IC(1,2)<0)
2084  IC(1,2)=(-1)*IC_12_corr;
2085  else
2086  IC(1,2)=IC_12_corr;
2087  }
2088  }
2089  }
2090  else{
2091  IC(1,2)=covAnyNumber_;
2092  }
2093 }
2094 //
2095 void CSCSegAlgoST::findDuplicates(std::vector<CSCSegment> & segments ){
2096  // this is intended for ME1/1a only - we have ghost segments because of the strips ganging
2097  // this function finds them (first the rechits by sharesInput() )
2098  // if a segment shares all the rechits with another segment it is a duplicate (even if
2099  // it has less rechits)
2100 
2101  for(std::vector<CSCSegment>::iterator it=segments.begin(); it != segments.end(); ++it) {
2102  std::vector<CSCSegment*> duplicateSegments;
2103  for(std::vector<CSCSegment>::iterator it2=segments.begin(); it2 != segments.end(); ++it2) {
2104  //
2105  bool allShared = true;
2106  if(it!=it2){
2107  allShared = it->sharesRecHits(*it2);
2108  }
2109  else{
2110  allShared = false;
2111  }
2112  //
2113  if(allShared){
2114  duplicateSegments.push_back(&(*it2));
2115  }
2116  }
2117  it->setDuplicateSegments(duplicateSegments);
2118  }
2119 
2120 }
2121 //
2122 
#define LogDebug(id)
double yweightPenaltyThreshold
Definition: CSCSegAlgoST.h:190
T getParameter(std::string const &) const
double protoChiUCorrection
Allow to correct the error matrix.
Definition: CSCSegAlgoST.h:199
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
float xx() const
Definition: LocalError.h:24
CSCSegment showerSeg(const CSCChamber *aChamber, ChamberHitContainer rechits)
bool preClustering_useChaining
Definition: CSCSegAlgoST.h:177
CSCDetId cscDetId() const
Definition: CSCRecHit2D.h:52
bool passCondNumber
The number to fource the Covariance.
Definition: CSCSegAlgoST.h:213
void doSlopesAndChi2(void)
std::vector< ChamberHitContainer > Psegments_noL5
Definition: CSCSegAlgoST.h:129
std::vector< float > curv_noL2_A
Definition: CSCSegAlgoST.h:143
double dXclusBoxMax
Definition: CSCSegAlgoST.h:173
void correctTheCovMatrix(CLHEP::HepMatrix &IC)
double_binary B
Definition: DDStreamer.cc:234
const CSCChamber * theChamber
Definition: CSCSegAlgoST.h:117
std::vector< float > weight_noL4_A
Definition: CSCSegAlgoST.h:137
void correctTheCovX(void)
std::vector< std::vector< const CSCRecHit2D * > > clusterHits(const CSCChamber *aChamber, ChamberHitContainer &rechits)
std::vector< float > curv_noL6_A
Definition: CSCSegAlgoST.h:147
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
void ChooseSegments2(int best_seg)
void ChooseSegments3(int best_seg)
std::vector< float > weight_noL3_B
Definition: CSCSegAlgoST.h:151
T y() const
Definition: PV3DBase.h:62
float protoSlope_v
Definition: CSCSegAlgoST.h:160
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:62
std::vector< ChamberHitContainer > Psegments_noL1
Definition: CSCSegAlgoST.h:125
CSCSegAlgoST(const edm::ParameterSet &ps)
Constructor.
Definition: CSCSegAlgoST.cc:33
bool isGoodToMerge(bool isME11a, ChamberHitContainer &newChain, ChamberHitContainer &oldChain)
std::vector< ChamberHitContainer > chosen_Psegments
Definition: CSCSegAlgoST.h:131
std::vector< ChamberHitContainer > Psegments
Definition: CSCSegAlgoST.h:123
ChamberHitContainer protoSegment
Definition: CSCSegAlgoST.h:158
int layer() const
Definition: CSCDetId.h:63
float protoSlope_u
Definition: CSCSegAlgoST.h:159
double double double z
double condSeed1_
Correct the error matrix for the condition number.
Definition: CSCSegAlgoST.h:209
LocalVector protoDirection
Definition: CSCSegAlgoST.h:164
std::vector< float > weight_noL3_A
Definition: CSCSegAlgoST.h:136
double hitDropLimit4Hits
Definition: CSCSegAlgoST.h:184
std::string chamberTypeName() const
float xy() const
Definition: LocalError.h:25
CSCSegAlgoShowering * showering_
Definition: CSCSegAlgoST.h:195
int nRecHits() const
Definition: CSCSegment.h:68
virtual ~CSCSegAlgoST()
Destructor.
Definition: CSCSegAlgoST.cc:91
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:42
CLHEP::HepMatrix AlgebraicMatrix
LocalError localPositionError() const
Definition: CSCRecHit2D.h:51
double BPMinImprovement
Definition: CSCSegAlgoST.h:180
float yy() const
Definition: LocalError.h:26
double covAnyNumber_
Allow to use any number for covariance for all RecHits.
Definition: CSCSegAlgoST.h:212
string unit
Definition: csvLumiCalc.py:46
T sqrt(T t)
Definition: SSEVec.h:46
std::vector< ChamberHitContainer > Psegments_noLx
Definition: CSCSegAlgoST.h:124
double curvePenaltyThreshold
Definition: CSCSegAlgoST.h:193
T z() const
Definition: PV3DBase.h:63
std::vector< ChamberHitContainer > Psegments_noL4
Definition: CSCSegAlgoST.h:128
tuple result
Definition: query.py:137
double protoNDF
Definition: CSCSegAlgoST.h:163
double hitDropLimit6Hits
Definition: CSCSegAlgoST.h:186
double protoChi2
Definition: CSCSegAlgoST.h:162
std::vector< float > weight_A
Definition: CSCSegAlgoST.h:132
std::vector< float > weight_noL1_B
Definition: CSCSegAlgoST.h:149
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:41
double chi2Norm_2D_
Definition: CSCSegAlgoST.h:201
float ChiSquaredProbability(double chiSquared, double nrDOF)
std::vector< CSCSegment > prune_bad_hits(const CSCChamber *aChamber, std::vector< CSCSegment > &segments)
void fitSlopes(void)
std::vector< float > curv_noL5_A
Definition: CSCSegAlgoST.h:146
std::vector< ChamberHitContainer > Psegments_noL6
Definition: CSCSegAlgoST.h:130
#define end
Definition: vmac.h:38
unsigned maxContrIndex
Chi^2 normalization for the initial fit.
Definition: CSCSegAlgoST.h:203
std::vector< const CSCRecHit2D * > ChamberHitContainer
Typedefs.
Definition: CSCSegAlgoST.h:40
bool onlyBestSegment
Definition: CSCSegAlgoST.h:181
void flipErrors(AlgebraicSymMatrix &) const
std::vector< ChamberHitContainer > Psegments_noL2
Definition: CSCSegAlgoST.h:126
std::vector< float > weight_noL5_A
Definition: CSCSegAlgoST.h:138
bool passCondNumber_2
Passage the condition number calculations.
Definition: CSCSegAlgoST.h:214
AlgebraicSymMatrix weightMatrix(void) const
std::vector< float > weight_noL1_A
Definition: CSCSegAlgoST.h:134
std::vector< float > curv_noL4_A
Definition: CSCSegAlgoST.h:145
std::vector< float > weight_noL2_A
Definition: CSCSegAlgoST.h:135
std::vector< float > chosen_weight_A
Definition: CSCSegAlgoST.h:140
Segments GoodSegments
Definition: CSCSegAlgoST.h:118
bool correctCov_
Correct the Error Matrix.
Definition: CSCSegAlgoST.h:198
int minHitsPerSegment
Definition: CSCSegAlgoST.h:170
ChamberHitContainer Psegments_hits
Definition: CSCSegAlgoST.h:121
double dYclusBoxMax
Definition: CSCSegAlgoST.h:174
void findDuplicates(std::vector< CSCSegment > &segments)
std::vector< float > curv_noL1_A
Definition: CSCSegAlgoST.h:142
void fillChiSquared(void)
AlgebraicSymMatrix calculateError(void) const
bool covToAnyNumberAll_
Allow to use any number for covariance (by hand)
Definition: CSCSegAlgoST.h:211
double b
Definition: hdecay.h:120
double theWeight(double coordinate_1, double coordinate_2, double coordinate_3, float layer_1, float layer_2, float layer_3)
Utility functions.
double curvePenalty
Definition: CSCSegAlgoST.h:194
std::vector< float > curv_A
Definition: CSCSegAlgoST.h:141
#define begin
Definition: vmac.h:31
bool prePrun_
The index of the worst x RecHit in Chi^2-X method.
Definition: CSCSegAlgoST.h:204
double chi2() const
Chi2 of the segment fit.
Definition: CSCSegment.h:58
double a
Definition: hdecay.h:121
std::vector< std::vector< const CSCRecHit2D * > > chainHits(const CSCChamber *aChamber, ChamberHitContainer &rechits)
bool preClustering
Definition: CSCSegAlgoST.h:176
double yweightPenalty
Definition: CSCSegAlgoST.h:191
CLHEP::HepSymMatrix AlgebraicSymMatrix
std::vector< float > weight_noL5_B
Definition: CSCSegAlgoST.h:153
std::vector< float > weight_B
Definition: CSCSegAlgoST.h:148
float a_yweightPenaltyThreshold[5][5]
Definition: CSCSegAlgoST.h:188
void ChooseSegments2a(std::vector< ChamberHitContainer > &best_segments, int best_seg)
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
std::vector< double > e_Cxx
Definition: CSCSegAlgoST.h:200
CLHEP::HepMatrix derivativeMatrix(void) const
void fillLocalDirection(void)
int maxRecHitsInCluster
Definition: CSCSegAlgoST.h:175
Definition: DDAxes.h:10
std::vector< CSCSegment > run(const CSCChamber *aChamber, ChamberHitContainer rechits)
Definition: CSCSegAlgoST.cc:96
std::vector< float > weight_noL6_A
Definition: CSCSegAlgoST.h:139
LocalPoint localPosition() const
Definition: CSCRecHit2D.h:50
std::vector< float > weight_noL2_B
Definition: CSCSegAlgoST.h:150
bool covToAnyNumber_
The correction parameters.
Definition: CSCSegAlgoST.h:210
double hitDropLimit5Hits
Definition: CSCSegAlgoST.h:185
std::vector< float > weight_noLx_A
Definition: CSCSegAlgoST.h:133
std::vector< float > curv_noL3_A
Definition: CSCSegAlgoST.h:144
std::vector< float > weight_noL6_B
Definition: CSCSegAlgoST.h:154
T x() const
Definition: PV3DBase.h:61
tuple size
Write out results.
mathSSE::Vec4< T > v
std::vector< CSCSegment > buildSegments(ChamberHitContainer rechits)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
ChamberHitContainer PAhits_onLayer[6]
Definition: CSCSegAlgoST.h:120
LocalPoint protoIntercept
Definition: CSCSegAlgoST.h:161
double condSeed2_
Definition: CSCSegAlgoST.h:209
std::vector< float > weight_noL4_B
Definition: CSCSegAlgoST.h:152
std::vector< ChamberHitContainer > Psegments_noL3
Definition: CSCSegAlgoST.h:127
double chi2Norm_3D_
Chi^2 normalization for the corrected fit.
Definition: CSCSegAlgoST.h:202
double prePrunLimit_
Definition: CSCSegAlgoST.h:207