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