00001
00010 #include "CSCSegAlgoST.h"
00011 #include "CSCSegAlgoShowering.h"
00012
00013 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00014
00015
00016
00017 #include "Geometry/CSCGeometry/interface/CSCLayer.h"
00018
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021
00022 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00023
00024 #include <algorithm>
00025 #include <cmath>
00026 #include <iostream>
00027 #include <string>
00028
00029
00030
00031
00032
00033 CSCSegAlgoST::CSCSegAlgoST(const edm::ParameterSet& ps) : CSCSegmentAlgorithm(ps), myName("CSCSegAlgoST") {
00034
00035 debug = ps.getUntrackedParameter<bool>("CSCDebug");
00036
00037
00038 minHitsPerSegment = ps.getParameter<int>("minHitsPerSegment");
00039
00040
00041 dXclusBoxMax = ps.getParameter<double>("dXclusBoxMax");
00042 dYclusBoxMax = ps.getParameter<double>("dYclusBoxMax");
00043 preClustering = ps.getParameter<bool>("preClustering");
00044 preClustering_useChaining = ps.getParameter<bool>("preClusteringUseChaining");
00045 Pruning = ps.getParameter<bool>("Pruning");
00046 BrutePruning = ps.getParameter<bool>("BrutePruning");
00047 BPMinImprovement = ps.getParameter<double>("BPMinImprovement");
00048
00049
00050
00051
00052 maxRecHitsInCluster = ps.getParameter<int>("maxRecHitsInCluster");
00053 onlyBestSegment = ps.getParameter<bool>("onlyBestSegment");
00054
00055 hitDropLimit4Hits = ps.getParameter<double>("hitDropLimit4Hits");
00056 hitDropLimit5Hits = ps.getParameter<double>("hitDropLimit5Hits");
00057 hitDropLimit6Hits = ps.getParameter<double>("hitDropLimit6Hits");
00058
00059 yweightPenaltyThreshold = ps.getParameter<double>("yweightPenaltyThreshold");
00060 yweightPenalty = ps.getParameter<double>("yweightPenalty");
00061
00062 curvePenaltyThreshold = ps.getParameter<double>("curvePenaltyThreshold");
00063 curvePenalty = ps.getParameter<double>("curvePenalty");
00064
00065 useShowering = ps.getParameter<bool>("useShowering");
00066 showering_ = new CSCSegAlgoShowering( ps );
00067
00069 correctCov_ = ps.getParameter<bool>("CorrectTheErrors");
00070 chi2Norm_2D_ = ps.getParameter<double>("NormChi2Cut2D");
00071 chi2Norm_3D_ = ps.getParameter<double>("NormChi2Cut3D");
00072 prePrun_ = ps.getParameter<bool>("prePrun");
00073 prePrunLimit_ = ps.getParameter<double>("prePrunLimit");
00074
00075 condSeed1_ = ps.getParameter<double>("SeedSmall");
00076 condSeed2_ = ps.getParameter<double>("SeedBig");
00077 covToAnyNumber_ = ps.getParameter<bool>("ForceCovariance");
00078 covToAnyNumberAll_ = ps.getParameter<bool>("ForceCovarianceAll");
00079 covAnyNumber_ = ps.getParameter<double>("Covariance");
00080 passCondNumber=false;
00081 passCondNumber_2=false;
00082 protoChiUCorrection=1.0;
00083 maxContrIndex=0;
00084 protoNDF = 1.;
00085
00086 }
00087
00088
00089
00090
00091 CSCSegAlgoST::~CSCSegAlgoST() {
00092 delete showering_;
00093 }
00094
00095
00096 std::vector<CSCSegment> CSCSegAlgoST::run(const CSCChamber* aChamber, ChamberHitContainer rechits) {
00097
00098
00099 theChamber = aChamber;
00100
00101 std::vector<CSCSegment> segments_temp;
00102 std::vector<CSCSegment> segments;
00103 std::vector<ChamberHitContainer> rechits_clusters;
00104
00105
00106
00107
00108 for(int a = 0; a<5; ++a) {
00109 for(int b = 0; b<5; ++b) {
00110 a_yweightPenaltyThreshold[a][b] = 0.0;
00111 }
00112 }
00113
00114 a_yweightPenaltyThreshold[1][1] = yweightPenaltyThreshold * 10.20;
00115 a_yweightPenaltyThreshold[1][2] = yweightPenaltyThreshold * 14.00;
00116 a_yweightPenaltyThreshold[1][3] = yweightPenaltyThreshold * 20.40;
00117 a_yweightPenaltyThreshold[1][4] = yweightPenaltyThreshold * 10.20;
00118 a_yweightPenaltyThreshold[2][1] = yweightPenaltyThreshold * 7.60;
00119 a_yweightPenaltyThreshold[2][2] = yweightPenaltyThreshold * 20.40;
00120 a_yweightPenaltyThreshold[3][1] = yweightPenaltyThreshold * 7.60;
00121 a_yweightPenaltyThreshold[3][2] = yweightPenaltyThreshold * 20.40;
00122 a_yweightPenaltyThreshold[4][1] = yweightPenaltyThreshold * 6.75;
00123
00124 if(preClustering) {
00125
00126 if(preClustering_useChaining){
00127
00128
00129
00130
00131
00132 rechits_clusters = chainHits( theChamber, rechits );
00133 }
00134 else{
00135
00136 rechits_clusters = clusterHits( theChamber, rechits );
00137 }
00138
00139 for(std::vector<ChamberHitContainer>::iterator sub_rechits = rechits_clusters.begin(); sub_rechits != rechits_clusters.end(); ++sub_rechits ) {
00140
00141 segments_temp.clear();
00142
00143 segments_temp = buildSegments( (*sub_rechits) );
00144
00145 segments.insert( segments.end(), segments_temp.begin(), segments_temp.end() );
00146 }
00147
00148 if( Pruning ) {
00149 segments_temp.clear();
00150 segments_temp = prune_bad_hits( theChamber, segments );
00151 segments.clear();
00152 segments.swap(segments_temp);
00153 }
00154 if ("ME1/a" == aChamber->specs()->chamberTypeName()){
00155 findDuplicates(segments);
00156 }
00157 return segments;
00158 }
00159 else {
00160 segments = buildSegments(rechits);
00161 if( Pruning ) {
00162 segments_temp.clear();
00163 segments_temp = prune_bad_hits( theChamber, segments );
00164 segments.clear();
00165 segments.swap(segments_temp);
00166 }
00167 if ("ME1/a" == aChamber->specs()->chamberTypeName()){
00168 findDuplicates(segments);
00169 }
00170 return segments;
00171
00172 }
00173 }
00174
00175
00176
00177
00178
00179
00180 std::vector<CSCSegment> CSCSegAlgoST::prune_bad_hits(const CSCChamber* aChamber, std::vector<CSCSegment> & segments) {
00181
00182
00183
00184
00185
00186 std::vector<CSCSegment> segments_temp;
00187 std::vector<ChamberHitContainer> rechits_clusters;
00188
00189 const float chi2ndfProbMin = 1.0e-4;
00190 bool use_brute_force = BrutePruning;
00191
00192 int hit_nr = 0;
00193 int hit_nr_worst = -1;
00194
00195
00196 for(std::vector<CSCSegment>::iterator it=segments.begin(); it != segments.end(); ++it) {
00197
00198
00199 if( (*it).nRecHits() <= minHitsPerSegment ) continue;
00200
00201 if( !use_brute_force ) {
00202
00203 float chisq = (*it).chi2();
00204 int nhits = (*it).nRecHits();
00205 LocalPoint localPos = (*it).localPosition();
00206 LocalVector segDir = (*it).localDirection();
00207 const CSCChamber* cscchamber = theChamber;
00208 float globZ ;
00209
00210 GlobalPoint globalPosition = cscchamber->toGlobal(localPos);
00211 globZ = globalPosition.z();
00212
00213
00214 if( ChiSquaredProbability((double)chisq,(double)(2*nhits-4)) < chi2ndfProbMin ) {
00215
00216
00217 std::vector<CSCRecHit2D> theseRecHits = (*it).specificRecHits();
00218 std::vector<CSCRecHit2D>::const_iterator iRH_worst;
00219
00220
00221 float xdist_local_worst_sig = -99999.;
00222 float xdist_local_2ndworst_sig = -99999.;
00223 float xdist_local_sig = -99999.;
00224
00225 hit_nr = 0;
00226 hit_nr_worst = -1;
00227
00228
00229 for ( std::vector<CSCRecHit2D>::const_iterator iRH = theseRecHits.begin(); iRH != theseRecHits.end(); ++iRH) {
00230
00231
00232
00233
00234 float loc_x_at_target ;
00235
00236
00237
00238
00239 loc_x_at_target = 0.;
00240
00241
00242
00243
00244
00245 const CSCLayer* csclayerRH = theChamber->layer((*iRH).cscDetId().layer());
00246 LocalPoint localPositionRH = (*iRH).localPosition();
00247 GlobalPoint globalPositionRH = csclayerRH->toGlobal(localPositionRH);
00248
00249 LocalError rerrlocal = (*iRH).localPositionError();
00250 float xxerr = rerrlocal.xx();
00251
00252 float target_z = globalPositionRH.z();
00253
00254 if(target_z > 0.) {
00255 loc_x_at_target = localPos.x() + (segDir.x()/fabs(segDir.z())*( target_z - globZ ));
00256
00257
00258 }
00259 else {
00260 loc_x_at_target = localPos.x() + ((-1)*segDir.x()/fabs(segDir.z())*( target_z - globZ ));
00261
00262
00263 }
00264
00265
00266
00267 xdist_local_sig = fabs((localPositionRH.x() -loc_x_at_target)/(xxerr));
00268
00269 if( xdist_local_sig > xdist_local_worst_sig ) {
00270 xdist_local_2ndworst_sig = xdist_local_worst_sig;
00271 xdist_local_worst_sig = xdist_local_sig;
00272 iRH_worst = iRH;
00273
00274 hit_nr_worst = hit_nr;
00275 }
00276 else if(xdist_local_sig > xdist_local_2ndworst_sig) {
00277 xdist_local_2ndworst_sig = xdist_local_sig;
00278
00279 }
00280 ++hit_nr;
00281 }
00282
00283
00284
00285
00286 if ( xdist_local_worst_sig / xdist_local_2ndworst_sig < 1.5 ) {
00287 hit_nr_worst = -1;
00288
00289 }
00290 }
00291 }
00292
00293
00294
00295
00296 std::vector< CSCRecHit2D > buffer;
00297 std::vector< std::vector< CSCRecHit2D > > reduced_segments;
00298 std::vector< CSCRecHit2D > theseRecHits = (*it).specificRecHits();
00299 float best_red_seg_prob = 0.0;
00300
00301 buffer.clear();
00302
00303 if( ChiSquaredProbability((double)(*it).chi2(),(double)((2*(*it).nRecHits())-4)) < chi2ndfProbMin ) {
00304
00305 buffer = theseRecHits;
00306
00307
00308
00309 if( use_brute_force ) {
00310 for(size_t bi = 0; bi < buffer.size(); ++bi) {
00311 reduced_segments.push_back(buffer);
00312 reduced_segments[bi].erase(reduced_segments[bi].begin()+(bi),reduced_segments[bi].begin()+(bi+1));
00313 }
00314 }
00315 else {
00316
00317 if( hit_nr_worst >= 0 && hit_nr_worst <= int(buffer.size()) ) {
00318
00319 buffer.erase(buffer.begin()+(hit_nr_worst),buffer.begin()+(hit_nr_worst+1));
00320 reduced_segments.push_back(buffer);
00321 }
00322 else {
00323
00324 reduced_segments.push_back(buffer);
00325 }
00326 }
00327 }
00328
00329
00330 for(size_t iSegment=0; iSegment<reduced_segments.size(); ++iSegment) {
00331
00332 protoSegment.clear();
00333 for(size_t m = 0; m<reduced_segments[iSegment].size(); ++m ) {
00334 protoSegment.push_back(&reduced_segments[iSegment][m]);
00335 }
00336 passCondNumber=false;
00337 passCondNumber_2 = false;
00338 protoChiUCorrection=1.0;
00339 doSlopesAndChi2();
00340
00341
00342
00343 if(correctCov_){
00344 if(protoChi2/protoNDF>chi2Norm_3D_){
00345 passCondNumber = true;
00346 doSlopesAndChi2();
00347 }
00348 if((protoChiUCorrection<1.00005)&&(protoChi2/protoNDF>chi2Norm_3D_)){
00349 passCondNumber_2=true;
00350 doSlopesAndChi2();
00351 }
00352 }
00353 fillLocalDirection();
00354
00355 AlgebraicSymMatrix protoErrors = calculateError();
00356
00357
00358 flipErrors( protoErrors );
00359
00360 CSCSegment temp(protoSegment, protoIntercept, protoDirection, protoErrors, protoChi2);
00361
00362
00363 if( ( ChiSquaredProbability((double)(*it).chi2(),(double)((2*(*it).nRecHits())-4))
00364 <
00365 (1./BPMinImprovement)*(ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4))) )
00366
00367 &&
00368 ( (ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4)))
00369 > best_red_seg_prob
00370 )
00371 &&
00372 ( (ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4))) > 1e-10 )
00373 ) {
00374 best_red_seg_prob = ChiSquaredProbability((double)temp.chi2(),(double)(2*temp.nRecHits()-4));
00375
00376
00377
00378 if( temp.nRecHits() >= minHitsPerSegment ) {
00379 (*it) = temp;
00380 }
00381 }
00382 }
00383 }
00384
00385 return segments;
00386
00387 }
00388
00389
00390
00391 std::vector< std::vector<const CSCRecHit2D*> > CSCSegAlgoST::clusterHits(const CSCChamber* aChamber, ChamberHitContainer & rechits) {
00392 theChamber = aChamber;
00393
00394 std::vector<ChamberHitContainer> rechits_clusters;
00395
00396
00397
00398
00399
00400 float dXclus_box = 0.0;
00401 float dYclus_box = 0.0;
00402
00403 std::vector<const CSCRecHit2D*> temp;
00404
00405 std::vector< ChamberHitContainer > seeds;
00406
00407 std::vector<float> running_meanX;
00408 std::vector<float> running_meanY;
00409
00410 std::vector<float> seed_minX;
00411 std::vector<float> seed_maxX;
00412 std::vector<float> seed_minY;
00413 std::vector<float> seed_maxY;
00414
00415
00416
00417
00418
00419
00420
00421
00422 for(unsigned int i = 0; i < rechits.size(); ++i) {
00423
00424 temp.clear();
00425
00426 temp.push_back(rechits[i]);
00427
00428 seeds.push_back(temp);
00429
00430
00431
00432
00433 running_meanX.push_back( rechits[i]->localPosition().x() );
00434 running_meanY.push_back( rechits[i]->localPosition().y() );
00435
00436
00437 seed_minX.push_back( rechits[i]->localPosition().x() );
00438 seed_maxX.push_back( rechits[i]->localPosition().x() );
00439 seed_minY.push_back( rechits[i]->localPosition().y() );
00440 seed_maxY.push_back( rechits[i]->localPosition().y() );
00441 }
00442
00443
00444
00445 for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
00446
00447 for(size_t MMM = NNN+1; MMM < seeds.size(); ++MMM) {
00448 if(running_meanX[MMM] == 999999. || running_meanX[NNN] == 999999. ) {
00449 LogDebug("CSCSegment|CSC") << "CSCSegmentST::clusterHits: Warning: Skipping used seeds, this should happen - inform developers!";
00450
00451 continue;
00452 }
00453
00454
00455
00456
00457
00458
00459 if ( running_meanX[NNN] > running_meanX[MMM] ) dXclus_box = seed_minX[NNN] - seed_maxX[MMM];
00460 else dXclus_box = seed_minX[MMM] - seed_maxX[NNN];
00461 if ( running_meanY[NNN] > running_meanY[MMM] ) dYclus_box = seed_minY[NNN] - seed_maxY[MMM];
00462 else dYclus_box = seed_minY[MMM] - seed_maxY[NNN];
00463
00464
00465 if( dXclus_box < dXclusBoxMax && dYclus_box < dYclusBoxMax ) {
00466
00467
00468
00469
00470 running_meanX[MMM] = (running_meanX[NNN]*seeds[NNN].size() + running_meanX[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
00471 running_meanY[MMM] = (running_meanY[NNN]*seeds[NNN].size() + running_meanY[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
00472
00473
00474 if ( seed_minX[NNN] <= seed_minX[MMM] ) seed_minX[MMM] = seed_minX[NNN];
00475 if ( seed_maxX[NNN] > seed_maxX[MMM] ) seed_maxX[MMM] = seed_maxX[NNN];
00476 if ( seed_minY[NNN] <= seed_minY[MMM] ) seed_minY[MMM] = seed_minY[NNN];
00477 if ( seed_maxY[NNN] > seed_maxY[MMM] ) seed_maxY[MMM] = seed_maxY[NNN];
00478
00479
00480 seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
00481
00482
00483 running_meanX[NNN] = 999999.;
00484 running_meanY[NNN] = 999999.;
00485
00486
00487 break;
00488 }
00489
00490 }
00491 }
00492
00493
00494
00495
00496 for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
00497 if(running_meanX[NNN] == 999999.) continue;
00498 rechits_clusters.push_back(seeds[NNN]);
00499 }
00500
00501
00502
00503 return rechits_clusters;
00504 }
00505
00506
00507 std::vector< std::vector<const CSCRecHit2D*> > CSCSegAlgoST::chainHits(const CSCChamber* aChamber, ChamberHitContainer & rechits) {
00508
00509 std::vector<ChamberHitContainer> rechits_chains;
00510
00511
00512 std::vector<const CSCRecHit2D*> temp;
00513
00514 std::vector< ChamberHitContainer > seeds;
00515
00516 std::vector <bool> usedCluster;
00517
00518
00519
00520
00521
00522 for(unsigned int i = 0; i < rechits.size(); ++i) {
00523 temp.clear();
00524 temp.push_back(rechits[i]);
00525 seeds.push_back(temp);
00526 usedCluster.push_back(false);
00527 }
00528 bool isME11a = false;
00529 if ("ME1/a" == aChamber->specs()->chamberTypeName()){
00530 isME11a = true;
00531 }
00532
00533 for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
00534 for(size_t MMM = NNN+1; MMM < seeds.size(); ++MMM) {
00535 if(usedCluster[MMM] || usedCluster[NNN]){
00536 continue;
00537 }
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549 bool goodToMerge = isGoodToMerge(isME11a, seeds[NNN], seeds[MMM]);
00550 if(goodToMerge){
00551
00552
00553
00554
00555 seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
00556
00557
00558 usedCluster[NNN] = true;
00559
00560
00561 break;
00562 }
00563
00564 }
00565 }
00566
00567
00568
00569
00570
00571 for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
00572 if(usedCluster[NNN]) continue;
00573 rechits_chains.push_back(seeds[NNN]);
00574 }
00575
00576
00577
00578 return rechits_chains;
00579 }
00580
00581 bool CSCSegAlgoST::isGoodToMerge(bool isME11a, ChamberHitContainer & newChain, ChamberHitContainer & oldChain) {
00582 for(size_t iRH_new = 0;iRH_new<newChain.size();++iRH_new){
00583 int layer_new = newChain[iRH_new]->cscDetId().layer()-1;
00584 int middleStrip_new = newChain[iRH_new]->nStrips()/2;
00585 int centralStrip_new = newChain[iRH_new]->channels(middleStrip_new);
00586 int centralWire_new = newChain[iRH_new]->hitWire();
00587 bool layerRequirementOK = false;
00588 bool stripRequirementOK = false;
00589 bool wireRequirementOK = false;
00590 bool goodToMerge = false;
00591 for(size_t iRH_old = 0;iRH_old<oldChain.size();++iRH_old){
00592 int layer_old = oldChain[iRH_old]->cscDetId().layer()-1;
00593 int middleStrip_old = oldChain[iRH_old]->nStrips()/2;
00594 int centralStrip_old = oldChain[iRH_old]->channels(middleStrip_old);
00595 int centralWire_old = oldChain[iRH_old]->hitWire();
00596
00597
00598
00599
00600
00601
00602
00603
00604 if(layer_new==layer_old+1 || layer_new==layer_old-1 ||
00605 layer_new==layer_old+2 || layer_new==layer_old-2 ||
00606 layer_new==layer_old+3 || layer_new==layer_old-3 ||
00607 layer_new==layer_old+4 || layer_new==layer_old-4 ){
00608 layerRequirementOK = true;
00609 }
00610 int allStrips = 48;
00611
00612
00613
00614 if(centralStrip_new==centralStrip_old ||
00615 centralStrip_new==centralStrip_old+1 || centralStrip_new==centralStrip_old-1 ||
00616 centralStrip_new==centralStrip_old+2 || centralStrip_new==centralStrip_old-2){
00617 stripRequirementOK = true;
00618 }
00619
00620 if(centralWire_new==centralWire_old ||
00621 centralWire_new==centralWire_old+1 || centralWire_new==centralWire_old-1 ||
00622 centralWire_new==centralWire_old+2 || centralWire_new==centralWire_old-2){
00623 wireRequirementOK = true;
00624 }
00625
00626 if(isME11a){
00627 if(centralStrip_new==centralStrip_old+1-allStrips || centralStrip_new==centralStrip_old-1-allStrips ||
00628 centralStrip_new==centralStrip_old+2-allStrips || centralStrip_new==centralStrip_old-2-allStrips ||
00629 centralStrip_new==centralStrip_old+1+allStrips || centralStrip_new==centralStrip_old-1+allStrips ||
00630 centralStrip_new==centralStrip_old+2+allStrips || centralStrip_new==centralStrip_old-2+allStrips){
00631 stripRequirementOK = true;
00632 }
00633 }
00634 if(layerRequirementOK && stripRequirementOK && wireRequirementOK){
00635 goodToMerge = true;
00636 return goodToMerge;
00637 }
00638 }
00639 }
00640 return false;
00641 }
00642
00643
00644
00645
00646 double CSCSegAlgoST::theWeight(double coordinate_1, double coordinate_2, double coordinate_3, float layer_1, float layer_2, float layer_3) {
00647 double sub_weight = 0;
00648 sub_weight = fabs(
00649 ( (coordinate_2 - coordinate_3) / (layer_2 - layer_3) ) -
00650 ( (coordinate_1 - coordinate_2) / (layer_1 - layer_2) )
00651 );
00652 return sub_weight;
00653 }
00654
00655
00656
00657
00658
00659 std::vector<CSCSegment> CSCSegAlgoST::buildSegments(ChamberHitContainer rechits) {
00660
00661
00662 std::vector<CSCSegment> segmentInChamber;
00663 segmentInChamber.clear();
00664
00665
00666 unsigned int thering = 999;
00667 unsigned int thestation = 999;
00668
00669
00670 std::vector<int> hits_onLayerNumber(6);
00671
00672 unsigned int UpperLimit = maxRecHitsInCluster;
00673 if (int(rechits.size()) < minHitsPerSegment) return segmentInChamber;
00674
00675 for(int iarray = 0; iarray <6; ++iarray) {
00676 PAhits_onLayer[iarray].clear();
00677 hits_onLayerNumber[iarray] = 0;
00678 }
00679
00680 chosen_Psegments.clear();
00681 chosen_weight_A.clear();
00682
00683 Psegments.clear();
00684 Psegments_noLx.clear();
00685 Psegments_noL1.clear();
00686 Psegments_noL2.clear();
00687 Psegments_noL3.clear();
00688 Psegments_noL4.clear();
00689 Psegments_noL5.clear();
00690 Psegments_noL6.clear();
00691
00692 Psegments_hits.clear();
00693
00694 weight_A.clear();
00695 weight_noLx_A.clear();
00696 weight_noL1_A.clear();
00697 weight_noL2_A.clear();
00698 weight_noL3_A.clear();
00699 weight_noL4_A.clear();
00700 weight_noL5_A.clear();
00701 weight_noL6_A.clear();
00702
00703 weight_B.clear();
00704 weight_noL1_B.clear();
00705 weight_noL2_B.clear();
00706 weight_noL3_B.clear();
00707 weight_noL4_B.clear();
00708 weight_noL5_B.clear();
00709 weight_noL6_B.clear();
00710
00711 curv_A.clear();
00712 curv_noL1_A.clear();
00713 curv_noL2_A.clear();
00714 curv_noL3_A.clear();
00715 curv_noL4_A.clear();
00716 curv_noL5_A.clear();
00717 curv_noL6_A.clear();
00718
00719
00720 int midlayer_pointer[6] = {0,0,2,3,3,4};
00721
00722
00723 int n_layers_occupied_tot = 0;
00724 int n_layers_processed = 0;
00725
00726 float min_weight_A = 99999.9;
00727 float min_weight_noLx_A = 99999.9;
00728
00729
00730
00731
00732
00733
00734
00735 int best_pseg = -1;
00736 int best_noLx_pseg = -1;
00737 int best_Layer_noLx = -1;
00738
00739
00740
00741
00742
00743
00744
00745
00746
00747
00748
00749
00750 for(size_t M = 0; M < rechits.size(); ++M) {
00751
00752 hits_onLayerNumber[ rechits[M]->cscDetId().layer()-1 ] += 1;
00753 if(hits_onLayerNumber[ rechits[M]->cscDetId().layer()-1 ] == 1 ) n_layers_occupied_tot += 1;
00754
00755 PAhits_onLayer[rechits[M]->cscDetId().layer()-1] .push_back(rechits[M]);
00756 }
00757
00758
00759
00760 int tothits = 0;
00761 int maxhits = 0;
00762 int nexthits = 0;
00763 int maxlayer = -1;
00764 int nextlayer = -1;
00765
00766 for(size_t i = 0; i< hits_onLayerNumber.size(); ++i){
00767
00768 tothits += hits_onLayerNumber[i];
00769 if (hits_onLayerNumber[i] > maxhits) {
00770 nextlayer = maxlayer;
00771 nexthits = maxhits;
00772 maxlayer = i;
00773 maxhits = hits_onLayerNumber[i];
00774 }
00775 else if (hits_onLayerNumber[i] > nexthits) {
00776 nextlayer = i;
00777 nexthits = hits_onLayerNumber[i];
00778 }
00779 }
00780
00781
00782 if (tothits > (int)UpperLimit) {
00783 if (n_layers_occupied_tot > 4) {
00784 tothits = tothits - hits_onLayerNumber[maxlayer];
00785 n_layers_occupied_tot = n_layers_occupied_tot - 1;
00786 PAhits_onLayer[maxlayer].clear();
00787 hits_onLayerNumber[maxlayer] = 0;
00788 }
00789 }
00790
00791 if (tothits > (int)UpperLimit) {
00792 if (n_layers_occupied_tot > 4) {
00793 tothits = tothits - hits_onLayerNumber[nextlayer];
00794 n_layers_occupied_tot = n_layers_occupied_tot - 1;
00795 PAhits_onLayer[nextlayer].clear();
00796 hits_onLayerNumber[nextlayer] = 0;
00797 }
00798 }
00799
00800 if (tothits > (int)UpperLimit){
00801
00802
00803
00804
00805 if (useShowering) {
00806 CSCSegment segShower = showering_->showerSeg(theChamber, rechits);
00807
00808
00809 if ( segShower.nRecHits() < 3 ) return segmentInChamber;
00810 if ( segShower.chi2() == -1 ) return segmentInChamber;
00811
00812 segmentInChamber.push_back(segShower);
00813 return segmentInChamber;
00814
00815 } else{
00816 LogDebug("CSCSegment|CSC") <<"Number of rechits in the cluster/chamber > "<< UpperLimit<<
00817 " ... Segment finding in the cluster/chamber canceled!";
00818
00819
00820 return segmentInChamber;
00821 }
00822 }
00823
00824
00825
00826
00827 if( rechits.size() > 0 ) {
00828 thering = rechits[0]->cscDetId().ring();
00829 thestation = rechits[0]->cscDetId().station();
00830
00831 }
00832
00833
00834
00835
00836 if( n_layers_occupied_tot < minHitsPerSegment ) {
00837 return segmentInChamber;
00838 }
00839
00840
00841
00842
00843
00844 for(int layer = 0; layer < 6; ++layer) {
00845
00846
00847
00848
00849
00850
00851
00852
00853 if( PAhits_onLayer[layer].size() > 0 ) {
00854 n_layers_processed += 1;
00855 }
00856
00857
00858 int orig_number_of_psegs = Psegments.size();
00859 int orig_number_of_noL1_psegs = Psegments_noL1.size();
00860 int orig_number_of_noL2_psegs = Psegments_noL2.size();
00861 int orig_number_of_noL3_psegs = Psegments_noL3.size();
00862 int orig_number_of_noL4_psegs = Psegments_noL4.size();
00863 int orig_number_of_noL5_psegs = Psegments_noL5.size();
00864 int orig_number_of_noL6_psegs = Psegments_noL6.size();
00865
00866
00867 for(int hit = 0; hit < int(PAhits_onLayer[layer].size()); ++hit) {
00868
00869
00870 if( orig_number_of_psegs == 0 ) {
00871
00872 Psegments_hits.push_back(PAhits_onLayer[layer][hit]);
00873
00874 Psegments.push_back(Psegments_hits);
00875 Psegments_noL6.push_back(Psegments_hits);
00876 Psegments_noL5.push_back(Psegments_hits);
00877 Psegments_noL4.push_back(Psegments_hits);
00878 Psegments_noL3.push_back(Psegments_hits);
00879 Psegments_noL2.push_back(Psegments_hits);
00880
00881
00882
00883 curv_A.push_back(0.0);
00884 curv_noL6_A.push_back(0.0);
00885 curv_noL5_A.push_back(0.0);
00886 curv_noL4_A.push_back(0.0);
00887 curv_noL3_A.push_back(0.0);
00888 curv_noL2_A.push_back(0.0);
00889
00890 weight_A.push_back(0.0);
00891 weight_noL6_A.push_back(0.0);
00892 weight_noL5_A.push_back(0.0);
00893 weight_noL4_A.push_back(0.0);
00894 weight_noL3_A.push_back(0.0);
00895 weight_noL2_A.push_back(0.0);
00896
00897 weight_B.push_back(0.0);
00898 weight_noL6_B.push_back(0.0);
00899 weight_noL5_B.push_back(0.0);
00900 weight_noL4_B.push_back(0.0);
00901 weight_noL3_B.push_back(0.0);
00902 weight_noL2_B.push_back(0.0);
00903
00904
00905 Psegments_hits .clear();
00906 }
00907 else {
00908 if( orig_number_of_noL1_psegs == 0 ) {
00909
00910 Psegments_hits.push_back(PAhits_onLayer[layer][hit]);
00911
00912 Psegments_noL1.push_back(Psegments_hits);
00913
00914
00915
00916 curv_noL1_A.push_back(0.0);
00917
00918 weight_noL1_A.push_back(0.0);
00919
00920 weight_noL1_B.push_back(0.0);
00921
00922
00923 Psegments_hits .clear();
00924
00925 }
00926
00927
00928
00929 for( int pseg = 0; pseg < orig_number_of_psegs; ++pseg ) {
00930
00931 int pseg_pos = (pseg)+((hit)*orig_number_of_psegs);
00932 int pseg_noL1_pos = (pseg)+((hit)*orig_number_of_noL1_psegs);
00933 int pseg_noL2_pos = (pseg)+((hit)*orig_number_of_noL2_psegs);
00934 int pseg_noL3_pos = (pseg)+((hit)*orig_number_of_noL3_psegs);
00935 int pseg_noL4_pos = (pseg)+((hit)*orig_number_of_noL4_psegs);
00936 int pseg_noL5_pos = (pseg)+((hit)*orig_number_of_noL5_psegs);
00937 int pseg_noL6_pos = (pseg)+((hit)*orig_number_of_noL6_psegs);
00938
00939
00940
00941
00942
00943 if( ! (hit == int(PAhits_onLayer[layer].size()-1)) ) {
00944
00945
00946 Psegments.push_back(Psegments[ pseg_pos ]);
00947 if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) Psegments_noL1.push_back(Psegments_noL1[ pseg_noL1_pos ]);
00948 if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) Psegments_noL2.push_back(Psegments_noL2[ pseg_noL2_pos ]);
00949 if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) Psegments_noL3.push_back(Psegments_noL3[ pseg_noL3_pos ]);
00950 if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) Psegments_noL4.push_back(Psegments_noL4[ pseg_noL4_pos ]);
00951 if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) Psegments_noL5.push_back(Psegments_noL5[ pseg_noL5_pos ]);
00952 if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) Psegments_noL6.push_back(Psegments_noL6[ pseg_noL6_pos ]);
00953
00954 weight_A.push_back(weight_A[ pseg_pos ]);
00955 if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) weight_noL1_A.push_back(weight_noL1_A[ pseg_noL1_pos ]);
00956 if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) weight_noL2_A.push_back(weight_noL2_A[ pseg_noL2_pos ]);
00957 if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) weight_noL3_A.push_back(weight_noL3_A[ pseg_noL3_pos ]);
00958 if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) weight_noL4_A.push_back(weight_noL4_A[ pseg_noL4_pos ]);
00959 if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) weight_noL5_A.push_back(weight_noL5_A[ pseg_noL5_pos ]);
00960 if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) weight_noL6_A.push_back(weight_noL6_A[ pseg_noL6_pos ]);
00961
00962 curv_A.push_back(curv_A[ pseg_pos ]);
00963 if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) curv_noL1_A.push_back(curv_noL1_A[ pseg_noL1_pos ]);
00964 if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) curv_noL2_A.push_back(curv_noL2_A[ pseg_noL2_pos ]);
00965 if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) curv_noL3_A.push_back(curv_noL3_A[ pseg_noL3_pos ]);
00966 if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) curv_noL4_A.push_back(curv_noL4_A[ pseg_noL4_pos ]);
00967 if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) curv_noL5_A.push_back(curv_noL5_A[ pseg_noL5_pos ]);
00968 if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) curv_noL6_A.push_back(curv_noL6_A[ pseg_noL6_pos ]);
00969
00970 weight_B.push_back(weight_B[ pseg_pos ]);
00971 if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) weight_noL1_B.push_back(weight_noL1_B[ pseg_noL1_pos ]);
00972 if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) weight_noL2_B.push_back(weight_noL2_B[ pseg_noL2_pos ]);
00973 if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) weight_noL3_B.push_back(weight_noL3_B[ pseg_noL3_pos ]);
00974 if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) weight_noL4_B.push_back(weight_noL4_B[ pseg_noL4_pos ]);
00975 if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) weight_noL5_B.push_back(weight_noL5_B[ pseg_noL5_pos ]);
00976 if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) weight_noL6_B.push_back(weight_noL6_B[ pseg_noL6_pos ]);
00977 }
00978
00979 Psegments[ pseg_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
00980 if (n_layers_processed != 2 && pseg < orig_number_of_noL1_psegs) Psegments_noL1[ pseg_noL1_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
00981 if (n_layers_processed != 2 && pseg < orig_number_of_noL2_psegs) Psegments_noL2[ pseg_noL2_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
00982 if (n_layers_processed != 3 && pseg < orig_number_of_noL3_psegs) Psegments_noL3[ pseg_noL3_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
00983 if (n_layers_processed != 4 && pseg < orig_number_of_noL4_psegs) Psegments_noL4[ pseg_noL4_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
00984 if (n_layers_processed != 5 && pseg < orig_number_of_noL5_psegs) Psegments_noL5[ pseg_noL5_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
00985 if (n_layers_processed != 6 && pseg < orig_number_of_noL6_psegs) Psegments_noL6[ pseg_noL6_pos ].push_back(PAhits_onLayer[ layer ][ hit ]);
00986
00987
00988
00989 if( Psegments[ pseg_pos ].size() > 2 ) {
00990
00991
00992
00993
00994 weight_A[ pseg_pos ] += theWeight(
00995 (*(Psegments[ pseg_pos ].end()-1 ))->localPosition().x(),
00996 (*(Psegments[ pseg_pos ].end()-2 ))->localPosition().x(),
00997 (*(Psegments[ pseg_pos ].end()-3 ))->localPosition().x(),
00998 float((*(Psegments[ pseg_pos ].end()-1))->cscDetId().layer()),
00999 float((*(Psegments[ pseg_pos ].end()-2))->cscDetId().layer()),
01000 float((*(Psegments[ pseg_pos ].end()-3))->cscDetId().layer())
01001 );
01002
01003 weight_B[ pseg_pos ] += theWeight(
01004 (*(Psegments[ pseg_pos ].end()-1 ))->localPosition().y(),
01005 (*(Psegments[ pseg_pos ].end()-2 ))->localPosition().y(),
01006 (*(Psegments[ pseg_pos ].end()-3 ))->localPosition().y(),
01007 float((*(Psegments[ pseg_pos ].end()-1))->cscDetId().layer()),
01008 float((*(Psegments[ pseg_pos ].end()-2))->cscDetId().layer()),
01009 float((*(Psegments[ pseg_pos ].end()-3))->cscDetId().layer())
01010 );
01011
01012
01013
01014 if(int(Psegments[ pseg_pos ].size()) == n_layers_occupied_tot) {
01015
01016 curv_A[ pseg_pos ] += theWeight(
01017 (*(Psegments[ pseg_pos ].end()-1 ))->localPosition().x(),
01018 (*(Psegments[ pseg_pos ].end()-midlayer_pointer[n_layers_occupied_tot-1] ))->localPosition().x(),
01019 (*(Psegments[ pseg_pos ].end()-n_layers_occupied_tot ))->localPosition().x(),
01020 float((*(Psegments[ pseg_pos ].end()-1))->cscDetId().layer()),
01021 float((*(Psegments[ pseg_pos ].end()-midlayer_pointer[n_layers_occupied_tot-1] ))->cscDetId().layer()),
01022 float((*(Psegments[ pseg_pos ].end()-n_layers_occupied_tot ))->cscDetId().layer())
01023 );
01024
01025 if (curv_A[ pseg_pos ] > curvePenaltyThreshold) weight_A[ pseg_pos ] = weight_A[ pseg_pos ] * curvePenalty;
01026
01027 if (weight_B[ pseg_pos ] > a_yweightPenaltyThreshold[thestation][thering]) weight_A[ pseg_pos ] = weight_A[ pseg_pos ] * yweightPenalty;
01028
01029 if (weight_A[ pseg_pos ] < min_weight_A ) {
01030 min_weight_A = weight_A[ pseg_pos ];
01031
01032
01033 best_pseg = pseg_pos ;
01034 }
01035
01036 }
01037
01038
01039
01040
01041 }
01042
01043 if ( n_layers_occupied_tot > 3 ) {
01044 if (pseg < orig_number_of_noL1_psegs && (n_layers_processed != 2)) {
01045 if(( Psegments_noL1[ pseg_noL1_pos ].size() > 2 ) ) {
01046
01047
01048
01049
01050 weight_noL1_A[ pseg_noL1_pos ] += theWeight(
01051 (*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->localPosition().x(),
01052 (*(Psegments_noL1[ pseg_noL1_pos ].end()-2 ))->localPosition().x(),
01053 (*(Psegments_noL1[ pseg_noL1_pos ].end()-3 ))->localPosition().x(),
01054 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-1))->cscDetId().layer()),
01055 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-2))->cscDetId().layer()),
01056 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-3))->cscDetId().layer())
01057 );
01058
01059 weight_noL1_B[ pseg_noL1_pos ] += theWeight(
01060 (*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->localPosition().y(),
01061 (*(Psegments_noL1[ pseg_noL1_pos ].end()-2 ))->localPosition().y(),
01062 (*(Psegments_noL1[ pseg_noL1_pos ].end()-3 ))->localPosition().y(),
01063 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-1))->cscDetId().layer()),
01064 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-2))->cscDetId().layer()),
01065 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-3))->cscDetId().layer())
01066 );
01067
01068
01069
01070 if(int(Psegments_noL1[ pseg_noL1_pos ].size()) == n_layers_occupied_tot -1 ) {
01071
01072 curv_noL1_A[ pseg_noL1_pos ] += theWeight(
01073 (*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->localPosition().x(),
01074 (*(Psegments_noL1[ pseg_noL1_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
01075 (*(Psegments_noL1[ pseg_noL1_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
01076 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-1 ))->cscDetId().layer()),
01077 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
01078 float((*(Psegments_noL1[ pseg_noL1_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
01079 );
01080
01081 if (curv_noL1_A[ pseg_noL1_pos ] > curvePenaltyThreshold) weight_noL1_A[ pseg_noL1_pos ] = weight_noL1_A[ pseg_noL1_pos ] * curvePenalty;
01082
01083 if (weight_noL1_B[ pseg_noL1_pos ] > a_yweightPenaltyThreshold[thestation][thering])
01084 weight_noL1_A[ pseg_noL1_pos ] = weight_noL1_A[ pseg_noL1_pos ] * yweightPenalty;
01085
01086 if (weight_noL1_A[ pseg_noL1_pos ] < min_weight_noLx_A ) {
01087 min_weight_noLx_A = weight_noL1_A[ pseg_noL1_pos ];
01088
01089
01090 best_noLx_pseg = pseg_noL1_pos;
01091 best_Layer_noLx = 1;
01092 }
01093
01094 }
01095
01096
01097
01098
01099 }
01100 }
01101 }
01102
01103 if ( n_layers_occupied_tot > 3 ) {
01104 if (pseg < orig_number_of_noL2_psegs && ( n_layers_processed != 2 )) {
01105 if(( Psegments_noL2[ pseg_noL2_pos ].size() > 2 )) {
01106
01107
01108
01109
01110 weight_noL2_A[ pseg_noL2_pos ] += theWeight(
01111 (*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->localPosition().x(),
01112 (*(Psegments_noL2[ pseg_noL2_pos ].end()-2 ))->localPosition().x(),
01113 (*(Psegments_noL2[ pseg_noL2_pos ].end()-3 ))->localPosition().x(),
01114 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-1))->cscDetId().layer()),
01115 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-2))->cscDetId().layer()),
01116 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-3))->cscDetId().layer())
01117 );
01118
01119 weight_noL2_B[ pseg_noL2_pos ] += theWeight(
01120 (*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->localPosition().y(),
01121 (*(Psegments_noL2[ pseg_noL2_pos ].end()-2 ))->localPosition().y(),
01122 (*(Psegments_noL2[ pseg_noL2_pos ].end()-3 ))->localPosition().y(),
01123 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-1))->cscDetId().layer()),
01124 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-2))->cscDetId().layer()),
01125 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-3))->cscDetId().layer())
01126 );
01127
01128
01129
01130 if(int(Psegments_noL2[ pseg_noL2_pos ].size()) == n_layers_occupied_tot -1 ) {
01131
01132 curv_noL2_A[ pseg_noL2_pos ] += theWeight(
01133 (*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->localPosition().x(),
01134 (*(Psegments_noL2[ pseg_noL2_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
01135 (*(Psegments_noL2[ pseg_noL2_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
01136 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-1 ))->cscDetId().layer()),
01137 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
01138 float((*(Psegments_noL2[ pseg_noL2_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
01139 );
01140
01141 if (curv_noL2_A[ pseg_noL2_pos ] > curvePenaltyThreshold) weight_noL2_A[ pseg_noL2_pos ] = weight_noL2_A[ pseg_noL2_pos ] * curvePenalty;
01142
01143 if (weight_noL2_B[ pseg_noL2_pos ] > a_yweightPenaltyThreshold[thestation][thering])
01144 weight_noL2_A[ pseg_noL2_pos ] = weight_noL2_A[ pseg_noL2_pos ] * yweightPenalty;
01145
01146 if (weight_noL2_A[ pseg_noL2_pos ] < min_weight_noLx_A ) {
01147 min_weight_noLx_A = weight_noL2_A[ pseg_noL2_pos ];
01148
01149
01150 best_noLx_pseg = pseg_noL2_pos;
01151 best_Layer_noLx = 2;
01152 }
01153
01154 }
01155
01156
01157
01158
01159 }
01160 }
01161 }
01162
01163 if ( n_layers_occupied_tot > 3 ) {
01164 if (pseg < orig_number_of_noL3_psegs && ( n_layers_processed != 3 )) {
01165 if(( Psegments_noL3[ pseg_noL3_pos ].size() > 2 )) {
01166
01167
01168
01169
01170 weight_noL3_A[ pseg_noL3_pos ] += theWeight(
01171 (*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->localPosition().x(),
01172 (*(Psegments_noL3[ pseg_noL3_pos ].end()-2 ))->localPosition().x(),
01173 (*(Psegments_noL3[ pseg_noL3_pos ].end()-3 ))->localPosition().x(),
01174 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-1))->cscDetId().layer()),
01175 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-2))->cscDetId().layer()),
01176 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-3))->cscDetId().layer())
01177 );
01178
01179 weight_noL3_B[ pseg_noL3_pos ] += theWeight(
01180 (*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->localPosition().y(),
01181 (*(Psegments_noL3[ pseg_noL3_pos ].end()-2 ))->localPosition().y(),
01182 (*(Psegments_noL3[ pseg_noL3_pos ].end()-3 ))->localPosition().y(),
01183 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-1))->cscDetId().layer()),
01184 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-2))->cscDetId().layer()),
01185 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-3))->cscDetId().layer())
01186 );
01187
01188
01189
01190 if(int(Psegments_noL3[ pseg_noL3_pos ].size()) == n_layers_occupied_tot -1 ) {
01191
01192 curv_noL3_A[ pseg_noL3_pos ] += theWeight(
01193 (*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->localPosition().x(),
01194 (*(Psegments_noL3[ pseg_noL3_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
01195 (*(Psegments_noL3[ pseg_noL3_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
01196 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-1 ))->cscDetId().layer()),
01197 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
01198 float((*(Psegments_noL3[ pseg_noL3_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
01199 );
01200
01201 if (curv_noL3_A[ pseg_noL3_pos ] > curvePenaltyThreshold) weight_noL3_A[ pseg_noL3_pos ] = weight_noL3_A[ pseg_noL3_pos ] * curvePenalty;
01202
01203 if (weight_noL3_B[ pseg_noL3_pos ] > a_yweightPenaltyThreshold[thestation][thering])
01204 weight_noL3_A[ pseg_noL3_pos ] = weight_noL3_A[ pseg_noL3_pos ] * yweightPenalty;
01205
01206 if (weight_noL3_A[ pseg_noL3_pos ] < min_weight_noLx_A ) {
01207 min_weight_noLx_A = weight_noL3_A[ pseg_noL3_pos ];
01208
01209
01210 best_noLx_pseg = pseg_noL3_pos;
01211 best_Layer_noLx = 3;
01212 }
01213
01214 }
01215
01216
01217
01218
01219 }
01220 }
01221 }
01222
01223 if ( n_layers_occupied_tot > 3 ) {
01224 if (pseg < orig_number_of_noL4_psegs && ( n_layers_processed != 4 )) {
01225 if(( Psegments_noL4[ pseg_noL4_pos ].size() > 2 )) {
01226
01227
01228
01229
01230 weight_noL4_A[ pseg_noL4_pos ] += theWeight(
01231 (*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->localPosition().x(),
01232 (*(Psegments_noL4[ pseg_noL4_pos ].end()-2 ))->localPosition().x(),
01233 (*(Psegments_noL4[ pseg_noL4_pos ].end()-3 ))->localPosition().x(),
01234 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-1))->cscDetId().layer()),
01235 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-2))->cscDetId().layer()),
01236 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-3))->cscDetId().layer())
01237 );
01238
01239 weight_noL4_B[ pseg_noL4_pos ] += theWeight(
01240 (*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->localPosition().y(),
01241 (*(Psegments_noL4[ pseg_noL4_pos ].end()-2 ))->localPosition().y(),
01242 (*(Psegments_noL4[ pseg_noL4_pos ].end()-3 ))->localPosition().y(),
01243 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-1))->cscDetId().layer()),
01244 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-2))->cscDetId().layer()),
01245 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-3))->cscDetId().layer())
01246 );
01247
01248
01249
01250 if(int(Psegments_noL4[ pseg_noL4_pos ].size()) == n_layers_occupied_tot -1 ) {
01251
01252 curv_noL4_A[ pseg_noL4_pos ] += theWeight(
01253 (*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->localPosition().x(),
01254 (*(Psegments_noL4[ pseg_noL4_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
01255 (*(Psegments_noL4[ pseg_noL4_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
01256 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-1 ))->cscDetId().layer()),
01257 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
01258 float((*(Psegments_noL4[ pseg_noL4_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
01259 );
01260
01261 if (curv_noL4_A[ pseg_noL4_pos ] > curvePenaltyThreshold) weight_noL4_A[ pseg_noL4_pos ] = weight_noL4_A[ pseg_noL4_pos ] * curvePenalty;
01262
01263 if (weight_noL4_B[ pseg_noL4_pos ] > a_yweightPenaltyThreshold[thestation][thering])
01264 weight_noL4_A[ pseg_noL4_pos ] = weight_noL4_A[ pseg_noL4_pos ] * yweightPenalty;
01265
01266 if (weight_noL4_A[ pseg_noL4_pos ] < min_weight_noLx_A ) {
01267 min_weight_noLx_A = weight_noL4_A[ pseg_noL4_pos ];
01268
01269
01270 best_noLx_pseg = pseg_noL4_pos;
01271 best_Layer_noLx = 4;
01272 }
01273
01274 }
01275
01276
01277
01278
01279 }
01280 }
01281 }
01282
01283 if ( n_layers_occupied_tot > 4 ) {
01284 if (pseg < orig_number_of_noL5_psegs && ( n_layers_processed != 5 )) {
01285 if(( Psegments_noL5[ pseg_noL5_pos ].size() > 2 )){
01286
01287
01288
01289
01290 weight_noL5_A[ pseg_noL5_pos ] += theWeight(
01291 (*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->localPosition().x(),
01292 (*(Psegments_noL5[ pseg_noL5_pos ].end()-2 ))->localPosition().x(),
01293 (*(Psegments_noL5[ pseg_noL5_pos ].end()-3 ))->localPosition().x(),
01294 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-1))->cscDetId().layer()),
01295 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-2))->cscDetId().layer()),
01296 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-3))->cscDetId().layer())
01297 );
01298
01299 weight_noL5_B[ pseg_noL5_pos ] += theWeight(
01300 (*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->localPosition().y(),
01301 (*(Psegments_noL5[ pseg_noL5_pos ].end()-2 ))->localPosition().y(),
01302 (*(Psegments_noL5[ pseg_noL5_pos ].end()-3 ))->localPosition().y(),
01303 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-1))->cscDetId().layer()),
01304 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-2))->cscDetId().layer()),
01305 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-3))->cscDetId().layer())
01306 );
01307
01308
01309
01310 if(int(Psegments_noL5[ pseg_noL5_pos ].size()) == n_layers_occupied_tot -1 ) {
01311
01312 curv_noL5_A[ pseg_noL5_pos ] += theWeight(
01313 (*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->localPosition().x(),
01314 (*(Psegments_noL5[ pseg_noL5_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
01315 (*(Psegments_noL5[ pseg_noL5_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
01316 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-1 ))->cscDetId().layer()),
01317 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
01318 float((*(Psegments_noL5[ pseg_noL5_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
01319 );
01320
01321 if (curv_noL5_A[ pseg_noL5_pos ] > curvePenaltyThreshold) weight_noL5_A[ pseg_noL5_pos ] = weight_noL5_A[ pseg_noL5_pos ] * curvePenalty;
01322
01323 if (weight_noL5_B[ pseg_noL5_pos ] > a_yweightPenaltyThreshold[thestation][thering])
01324 weight_noL5_A[ pseg_noL5_pos ] = weight_noL5_A[ pseg_noL5_pos ] * yweightPenalty;
01325
01326 if (weight_noL5_A[ pseg_noL5_pos ] < min_weight_noLx_A ) {
01327 min_weight_noLx_A = weight_noL5_A[ pseg_noL5_pos ];
01328
01329
01330 best_noLx_pseg = pseg_noL5_pos;
01331 best_Layer_noLx = 5;
01332 }
01333
01334 }
01335
01336
01337
01338
01339 }
01340 }
01341 }
01342
01343 if ( n_layers_occupied_tot > 5 ) {
01344 if (pseg < orig_number_of_noL6_psegs && ( n_layers_processed != 6 )) {
01345 if(( Psegments_noL6[ pseg_noL6_pos ].size() > 2 )){
01346
01347
01348
01349
01350 weight_noL6_A[ pseg_noL6_pos ] += theWeight(
01351 (*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->localPosition().x(),
01352 (*(Psegments_noL6[ pseg_noL6_pos ].end()-2 ))->localPosition().x(),
01353 (*(Psegments_noL6[ pseg_noL6_pos ].end()-3 ))->localPosition().x(),
01354 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-1))->cscDetId().layer()),
01355 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-2))->cscDetId().layer()),
01356 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-3))->cscDetId().layer())
01357 );
01358
01359 weight_noL6_B[ pseg_noL6_pos ] += theWeight(
01360 (*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->localPosition().y(),
01361 (*(Psegments_noL6[ pseg_noL6_pos ].end()-2 ))->localPosition().y(),
01362 (*(Psegments_noL6[ pseg_noL6_pos ].end()-3 ))->localPosition().y(),
01363 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-1))->cscDetId().layer()),
01364 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-2))->cscDetId().layer()),
01365 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-3))->cscDetId().layer())
01366 );
01367
01368
01369
01370 if(int(Psegments_noL6[ pseg_noL6_pos ].size()) == n_layers_occupied_tot -1 ) {
01371
01372 curv_noL6_A[ pseg_noL6_pos ] += theWeight(
01373 (*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->localPosition().x(),
01374 (*(Psegments_noL6[ pseg_noL6_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->localPosition().x(),
01375 (*(Psegments_noL6[ pseg_noL6_pos ].end()-(n_layers_occupied_tot-1) ))->localPosition().x(),
01376 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-1 ))->cscDetId().layer()),
01377 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-midlayer_pointer[n_layers_occupied_tot-2] ))->cscDetId().layer()),
01378 float((*(Psegments_noL6[ pseg_noL6_pos ].end()-(n_layers_occupied_tot-1) ))->cscDetId().layer())
01379 );
01380
01381 if (curv_noL6_A[ pseg_noL6_pos ] > curvePenaltyThreshold) weight_noL6_A[ pseg_noL6_pos ] = weight_noL6_A[ pseg_noL6_pos ] * curvePenalty;
01382
01383 if (weight_noL6_B[ pseg_noL6_pos ] > a_yweightPenaltyThreshold[thestation][thering])
01384 weight_noL6_A[ pseg_noL6_pos ] = weight_noL6_A[ pseg_noL6_pos ] * yweightPenalty;
01385
01386 if (weight_noL6_A[ pseg_noL6_pos ] < min_weight_noLx_A ) {
01387 min_weight_noLx_A = weight_noL6_A[ pseg_noL6_pos ];
01388
01389
01390 best_noLx_pseg = pseg_noL6_pos;
01391 best_Layer_noLx = 6;
01392 }
01393
01394 }
01395
01396
01397
01398
01399 }
01400 }
01401 }
01402
01403 }
01404 }
01405 }
01406 }
01407
01408
01409
01410
01411
01412
01413
01414
01415
01416
01417
01418
01419
01420
01421
01422
01423
01424
01425
01426
01427 int chosen_pseg = best_pseg;
01428 if (best_pseg<0) {
01429 return segmentInChamber;
01430 }
01431 chosen_Psegments = (Psegments);
01432 chosen_weight_A = (weight_A);
01433
01434 float hit_drop_limit = -999999.999;
01435
01436
01437 switch ( n_layers_processed ) {
01438 case 1 :
01439
01440 break;
01441 case 2 :
01442
01443 break;
01444 case 3 :
01445
01446 break;
01447 case 4 :
01448 hit_drop_limit = hitDropLimit6Hits * (1./2.) * hitDropLimit4Hits;
01449 if ((best_Layer_noLx < 1) || (best_Layer_noLx > 4)) {
01450
01451 }
01452 if ((best_Layer_noLx == 2) || (best_Layer_noLx == 3)) hit_drop_limit = hit_drop_limit * (1./2.);
01453 break;
01454 case 5 :
01455 hit_drop_limit = hitDropLimit6Hits * (2./3.) * hitDropLimit5Hits;
01456 if ((best_Layer_noLx < 1) || (best_Layer_noLx > 5)) {
01457
01458 }
01459 if ((best_Layer_noLx == 2) || (best_Layer_noLx == 4)) hit_drop_limit = hit_drop_limit * (1./2.);
01460 if (best_Layer_noLx == 3) hit_drop_limit = hit_drop_limit * (1./3.);
01461 break;
01462 case 6 :
01463 hit_drop_limit = hitDropLimit6Hits * (3./4.);
01464 if ((best_Layer_noLx < 1) || (best_Layer_noLx > 6)) {
01465
01466 }
01467 if ((best_Layer_noLx == 2) || (best_Layer_noLx == 5)) hit_drop_limit = hit_drop_limit * (1./2.);
01468 if ((best_Layer_noLx == 3) || (best_Layer_noLx == 4)) hit_drop_limit = hit_drop_limit * (1./3.);
01469 break;
01470
01471 default :
01472
01473 LogDebug("CSCSegment|CSC") <<"CSCSegAlgoST: Unexpected number of layers with hits - please inform developers.";
01474
01475 hit_drop_limit = 0.1;
01476 }
01477
01478
01479 switch ( best_Layer_noLx ) {
01480 case 1 :
01481 Psegments_noLx.clear();
01482 Psegments_noLx = Psegments_noL1;
01483 weight_noLx_A.clear();
01484 weight_noLx_A = weight_noL1_A;
01485 break;
01486 case 2 :
01487 Psegments_noLx.clear();
01488 Psegments_noLx = Psegments_noL2;
01489 weight_noLx_A.clear();
01490 weight_noLx_A = weight_noL2_A;
01491 break;
01492 case 3 :
01493 Psegments_noLx.clear();
01494 Psegments_noLx = Psegments_noL3;
01495 weight_noLx_A.clear();
01496 weight_noLx_A = weight_noL3_A;
01497 break;
01498 case 4 :
01499 Psegments_noLx.clear();
01500 Psegments_noLx = Psegments_noL4;
01501 weight_noLx_A.clear();
01502 weight_noLx_A = weight_noL4_A;
01503 break;
01504 case 5 :
01505 Psegments_noLx.clear();
01506 Psegments_noLx = Psegments_noL5;
01507 weight_noLx_A.clear();
01508 weight_noLx_A = weight_noL5_A;
01509 break;
01510 case 6 :
01511 Psegments_noLx.clear();
01512 Psegments_noLx = Psegments_noL6;
01513 weight_noLx_A.clear();
01514 weight_noLx_A = weight_noL6_A;
01515 break;
01516
01517 default :
01518
01519 Psegments_noLx.clear();
01520 weight_noLx_A.clear();
01521 }
01522
01523 if( min_weight_A > 0. ) {
01524 if ( min_weight_noLx_A/min_weight_A < hit_drop_limit ) {
01525
01526
01527
01528
01529 chosen_pseg = best_noLx_pseg;
01530 chosen_Psegments.clear();
01531 chosen_weight_A.clear();
01532 chosen_Psegments = (Psegments_noLx);
01533 chosen_weight_A = (weight_noLx_A);
01534 }
01535 }
01536
01537 if(onlyBestSegment) {
01538 ChooseSegments2a( chosen_Psegments, chosen_pseg );
01539 }
01540 else {
01541 ChooseSegments3( chosen_Psegments, chosen_weight_A, chosen_pseg );
01542 }
01543
01544 for(unsigned int iSegment=0; iSegment<GoodSegments.size();++iSegment){
01545 protoSegment = GoodSegments[iSegment];
01546 passCondNumber=false;
01547 passCondNumber_2 = false;
01548 protoChiUCorrection=1.0;
01549 doSlopesAndChi2();
01550
01551
01552
01553
01554 if(correctCov_){
01555
01556
01557
01558
01559
01560 if(protoChi2/protoNDF>chi2Norm_3D_){
01561 passCondNumber = true;
01562 doSlopesAndChi2();
01563 }
01564 if(protoChiUCorrection<1.00005){
01565 LogDebug("CSCSegment|segmWierd") << "Wierd segment, ErrXX scaled, refit " <<std::endl;
01566 if(protoChi2/protoNDF>chi2Norm_3D_){
01567
01568
01569
01570
01571
01572 LogDebug("CSCSegment|segmWierd") << "Wierd segment, ErrXY changed to match cond. number, refit " << std::endl;
01573 passCondNumber_2=true;
01574 doSlopesAndChi2();
01575 }
01576 }
01577
01578
01579
01580
01581
01582 if(prePrun_ && (sqrt(protoChiUCorrection)>prePrunLimit_) &&
01583 (protoSegment.size()>3)){
01584 LogDebug("CSCSegment|segmWierd") << "Scale factor protoChiUCorrection too big, pre-Prune, refit " << std::endl;
01585 protoSegment.erase(protoSegment.begin()+(maxContrIndex),
01586 protoSegment.begin()+(maxContrIndex+1));
01587 doSlopesAndChi2();
01588 }
01589 }
01590
01591 fillLocalDirection();
01592
01593 AlgebraicSymMatrix protoErrors = calculateError();
01594
01595
01596 flipErrors( protoErrors );
01597
01598 CSCSegment temp(protoSegment, protoIntercept, protoDirection, protoErrors, protoChi2);
01599 segmentInChamber.push_back(temp);
01600 }
01601 return segmentInChamber;
01602 }
01603
01604 void CSCSegAlgoST::ChooseSegments2a(std::vector< ChamberHitContainer > & chosen_segments, int chosen_seg) {
01605
01606 GoodSegments.clear();
01607 GoodSegments.push_back( chosen_segments[chosen_seg] );
01608 }
01609
01610 void CSCSegAlgoST::ChooseSegments3(std::vector< ChamberHitContainer > & chosen_segments, std::vector< float > & chosen_weight, int chosen_seg) {
01611
01612 int SumCommonHits = 0;
01613 GoodSegments.clear();
01614 int nr_remaining_candidates;
01615 unsigned int nr_of_segment_candidates;
01616
01617 nr_remaining_candidates = nr_of_segment_candidates = chosen_segments.size();
01618
01619
01620 GoodSegments.push_back( chosen_segments[ chosen_seg ] );
01621
01622 float chosen_weight_temp = 999999.;
01623 int chosen_seg_temp = -1;
01624
01625
01626 while( nr_remaining_candidates > 0 ) {
01627
01628 for(unsigned int iCand=0; iCand < nr_of_segment_candidates; ++iCand) {
01629
01630 if( chosen_weight[iCand] < 0. ) continue;
01631 SumCommonHits = 0;
01632
01633 for( int ihits = 0; ihits < int(chosen_segments[iCand].size()); ++ihits ) {
01634 if( chosen_segments[iCand][ihits] == chosen_segments[chosen_seg][ihits]) {
01635 ++SumCommonHits;
01636 }
01637 }
01638
01639
01640 if(SumCommonHits>1) {
01641 chosen_weight[iCand] = -1.;
01642 nr_remaining_candidates -= 1;
01643 }
01644 else {
01645
01646 if( chosen_weight[ iCand ] < chosen_weight_temp ) {
01647 chosen_weight_temp = chosen_weight[ iCand ];
01648 chosen_seg_temp = iCand ;
01649 }
01650 }
01651 }
01652
01653 if( chosen_seg_temp > -1 ) GoodSegments.push_back( chosen_segments[ chosen_seg_temp ] );
01654
01655 chosen_seg = chosen_seg_temp;
01656
01657 chosen_weight_temp = 999999;
01658 chosen_seg_temp = -1;
01659 }
01660 }
01661
01662 void CSCSegAlgoST::ChooseSegments2(int best_seg) {
01663
01664 std::vector <unsigned int> BadCandidate;
01665 int SumCommonHits =0;
01666 GoodSegments.clear();
01667 BadCandidate.clear();
01668 for(unsigned int iCand=0;iCand<Psegments.size();++iCand) {
01669
01670 for(unsigned int iiCand=iCand+1;iiCand<Psegments.size();++iiCand){
01671
01672 SumCommonHits =0;
01673 if( Psegments[iCand].size() != Psegments[iiCand].size() ) {
01674 LogDebug("CSCSegment|CSC") <<"CSCSegmentST::ChooseSegments2: ALARM!! THIS should not happen!!";
01675
01676 }
01677 else {
01678 for( int ihits = 0; ihits < int(Psegments[iCand].size()); ++ihits ) {
01679 if( Psegments[iCand][ihits] == Psegments[iiCand][ihits]) {
01680 ++SumCommonHits;
01681 }
01682 }
01683 }
01684 if(SumCommonHits>1) {
01685 if( weight_A[iCand]>weight_A[iiCand] ) {
01686 BadCandidate.push_back(iCand);
01687
01688 }
01689 else{
01690 BadCandidate.push_back(iiCand);
01691
01692 }
01693 }
01694 }
01695 }
01696 bool discard;
01697 for(unsigned int isegm=0;isegm<Psegments.size();++isegm) {
01698
01699
01700 discard = false;
01701 for(unsigned int ibad=0;ibad<BadCandidate.size();++ibad) {
01702
01703 if(isegm == BadCandidate[ibad]) {
01704 discard = true;
01705 }
01706 }
01707 if(!discard) {
01708 GoodSegments.push_back( Psegments[isegm] );
01709 }
01710 }
01711 }
01712
01713
01714
01715
01716 void CSCSegAlgoST::doSlopesAndChi2(){
01717 fitSlopes();
01718 fillChiSquared();
01719 }
01720
01721
01722
01723
01724
01725 void CSCSegAlgoST::fitSlopes() {
01726 e_Cxx.clear();
01727 if(passCondNumber && !passCondNumber_2){
01728 correctTheCovX();
01729 if(e_Cxx.size()!=protoSegment.size()){
01730 LogDebug("CSCSegment|segmWierd") << "e_Cxx.size()!=protoSegment.size() IT IS A SERIOUS PROBLEM!!! " <<std::endl;
01731 }
01732 }
01733 CLHEP::HepMatrix M(4,4,0);
01734 CLHEP::HepVector B(4,0);
01735 ChamberHitContainer::const_iterator ih = protoSegment.begin();
01736 for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
01737 const CSCRecHit2D& hit = (**ih);
01738 const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
01739 GlobalPoint gp = layer->toGlobal(hit.localPosition());
01740 LocalPoint lp = theChamber->toLocal(gp);
01741
01742 double u = lp.x();
01743 double v = lp.y();
01744 double z = lp.z();
01745
01746 CLHEP::HepMatrix IC(2,2);
01747 if(passCondNumber&& !passCondNumber_2){
01748 IC(1,1) = e_Cxx.at(ih-protoSegment.begin());
01749 }
01750 else{
01751 IC(1,1) = hit.localPositionError().xx();
01752 }
01753
01754 IC(1,2) = hit.localPositionError().xy();
01755 IC(2,2) = hit.localPositionError().yy();
01756 IC(2,1) = IC(1,2);
01758 if(passCondNumber_2){
01759 correctTheCovMatrix(IC);
01760 }
01761
01762 int ierr = 0;
01763 IC.invert(ierr);
01764 if (ierr != 0) {
01765 LogDebug("CSCSegment|CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC;
01766
01767 }
01768
01769 M(1,1) += IC(1,1);
01770 M(1,2) += IC(1,2);
01771 M(1,3) += IC(1,1) * z;
01772 M(1,4) += IC(1,2) * z;
01773 B(1) += u * IC(1,1) + v * IC(1,2);
01774
01775 M(2,1) += IC(2,1);
01776 M(2,2) += IC(2,2);
01777 M(2,3) += IC(2,1) * z;
01778 M(2,4) += IC(2,2) * z;
01779 B(2) += u * IC(2,1) + v * IC(2,2);
01780
01781 M(3,1) += IC(1,1) * z;
01782 M(3,2) += IC(1,2) * z;
01783 M(3,3) += IC(1,1) * z * z;
01784 M(3,4) += IC(1,2) * z * z;
01785 B(3) += ( u * IC(1,1) + v * IC(1,2) ) * z;
01786
01787 M(4,1) += IC(2,1) * z;
01788 M(4,2) += IC(2,2) * z;
01789 M(4,3) += IC(2,1) * z * z;
01790 M(4,4) += IC(2,2) * z * z;
01791 B(4) += ( u * IC(2,1) + v * IC(2,2) ) * z;
01792 }
01793 CLHEP::HepVector p = solve(M, B);
01794
01795
01796
01797 protoIntercept = LocalPoint(p(1), p(2), 0.);
01798 protoSlope_u = p(3);
01799 protoSlope_v = p(4);
01800 }
01801
01802
01803
01804
01805
01806 void CSCSegAlgoST::fillChiSquared() {
01807
01808 double chsq = 0.;
01809
01810 ChamberHitContainer::const_iterator ih;
01811 for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
01812
01813 const CSCRecHit2D& hit = (**ih);
01814 const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
01815 GlobalPoint gp = layer->toGlobal(hit.localPosition());
01816 LocalPoint lp = theChamber->toLocal(gp);
01817
01818 double u = lp.x();
01819 double v = lp.y();
01820 double z = lp.z();
01821
01822 double du = protoIntercept.x() + protoSlope_u * z - u;
01823 double dv = protoIntercept.y() + protoSlope_v * z - v;
01824
01825 CLHEP::HepMatrix IC(2,2);
01826 if(passCondNumber&& !passCondNumber_2){
01827 IC(1,1) = e_Cxx.at(ih-protoSegment.begin());
01828 }
01829 else{
01830 IC(1,1) = hit.localPositionError().xx();
01831 }
01832
01833 IC(1,2) = hit.localPositionError().xy();
01834 IC(2,2) = hit.localPositionError().yy();
01835 IC(2,1) = IC(1,2);
01837 if(passCondNumber_2){
01838 correctTheCovMatrix(IC);
01839 }
01840
01841
01842 int ierr = 0;
01843 IC.invert(ierr);
01844 if (ierr != 0) {
01845 LogDebug("CSCSegment|CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC;
01846
01847
01848 }
01849
01850 chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
01851 }
01852
01853 protoChi2 = chsq;
01854 protoNDF = 2.*protoSegment.size() - 4;
01855 }
01856
01857
01858
01859 void CSCSegAlgoST::fillLocalDirection() {
01860
01861
01862
01863 double dxdz = protoSlope_u;
01864 double dydz = protoSlope_v;
01865 double dz = 1./sqrt(1. + dxdz*dxdz + dydz*dydz);
01866 double dx = dz*dxdz;
01867 double dy = dz*dydz;
01868 LocalVector localDir(dx,dy,dz);
01869
01870
01871
01872
01873
01874 double globalZpos = ( theChamber->toGlobal( protoIntercept ) ).z();
01875 double globalZdir = ( theChamber->toGlobal( localDir ) ).z();
01876 double directionSign = globalZpos * globalZdir;
01877 protoDirection = (directionSign * localDir).unit();
01878 }
01879
01880
01881
01882 AlgebraicSymMatrix CSCSegAlgoST::weightMatrix() const {
01883
01884 std::vector<const CSCRecHit2D*>::const_iterator it;
01885 int nhits = protoSegment.size();
01886 AlgebraicSymMatrix matrix(2*nhits, 0);
01887 int row = 0;
01888
01889 for (it = protoSegment.begin(); it != protoSegment.end(); ++it) {
01890
01891 const CSCRecHit2D& hit = (**it);
01892 ++row;
01893 matrix(row, row) = protoChiUCorrection*hit.localPositionError().xx();
01894 matrix(row, row+1) = hit.localPositionError().xy();
01895 ++row;
01896 matrix(row, row-1) = hit.localPositionError().xy();
01897 matrix(row, row) = hit.localPositionError().yy();
01898 }
01899 int ierr;
01900 matrix.invert(ierr);
01901 return matrix;
01902 }
01903
01904
01905
01906
01907
01908 CLHEP::HepMatrix CSCSegAlgoST::derivativeMatrix() const {
01909
01910 ChamberHitContainer::const_iterator it;
01911 int nhits = protoSegment.size();
01912 CLHEP::HepMatrix matrix(2*nhits, 4);
01913 int row = 0;
01914
01915 for(it = protoSegment.begin(); it != protoSegment.end(); ++it) {
01916
01917 const CSCRecHit2D& hit = (**it);
01918 const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
01919 GlobalPoint gp = layer->toGlobal(hit.localPosition());
01920 LocalPoint lp = theChamber->toLocal(gp);
01921 float z = lp.z();
01922 ++row;
01923 matrix(row, 1) = 1.;
01924 matrix(row, 3) = z;
01925 ++row;
01926 matrix(row, 2) = 1.;
01927 matrix(row, 4) = z;
01928 }
01929 return matrix;
01930 }
01931
01932
01933
01934
01935
01936 AlgebraicSymMatrix CSCSegAlgoST::calculateError() const {
01937
01938 AlgebraicSymMatrix weights = weightMatrix();
01939 AlgebraicMatrix A = derivativeMatrix();
01940
01941
01942
01943 int ierr;
01944 AlgebraicSymMatrix result = weights.similarityT(A);
01945 result.invert(ierr);
01946
01947
01948 return result;
01949 }
01950
01951
01952 void CSCSegAlgoST::flipErrors( AlgebraicSymMatrix& a ) const {
01953
01954
01955
01956
01957 AlgebraicSymMatrix hold( a );
01958
01959
01960 a(1,1) = hold(3,3);
01961 a(1,2) = hold(3,4);
01962 a(2,1) = hold(4,3);
01963 a(2,2) = hold(4,4);
01964
01965
01966 a(3,3) = hold(1,1);
01967 a(3,4) = hold(1,2);
01968 a(4,3) = hold(2,1);
01969 a(4,4) = hold(2,2);
01970
01971
01972 a(4,1) = hold(2,3);
01973 a(3,2) = hold(1,4);
01974 a(2,3) = hold(4,1);
01975 a(1,4) = hold(3,2);
01976 }
01977
01978 void CSCSegAlgoST::correctTheCovX(void){
01979 std::vector<double> uu, vv, zz;
01980
01981 e_Cxx.clear();
01982 double sum_U_err=0.0;
01983 double sum_Z_U_err=0.0;
01984 double sum_Z2_U_err=0.0;
01985 double sum_U_U_err=0.0;
01986 double sum_UZ_U_err=0.0;
01987 std::vector<double> chiUZind;
01988 std::vector<double>::iterator chiContribution;
01989 double chiUZ=0.0;
01990 ChamberHitContainer::const_iterator ih = protoSegment.begin();
01991 for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
01992 const CSCRecHit2D& hit = (**ih);
01993 e_Cxx.push_back(hit.localPositionError().xx());
01994
01995 const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
01996 GlobalPoint gp = layer->toGlobal(hit.localPosition());
01997 LocalPoint lp = theChamber->toLocal(gp);
01998
01999 double u = lp.x();
02000 double v = lp.y();
02001 double z = lp.z();
02002 uu.push_back(u);
02003 vv.push_back(v);
02004 zz.push_back(z);
02006 sum_U_err += 1./e_Cxx.back();
02007 sum_Z_U_err += z/e_Cxx.back();
02008 sum_Z2_U_err += (z*z)/e_Cxx.back();
02009 sum_U_U_err += u/e_Cxx.back();
02010 sum_UZ_U_err += (u*z)/e_Cxx.back();
02011 }
02012
02015
02016 double denom=sum_U_err*sum_Z2_U_err-pow(sum_Z_U_err,2);
02017 double U0=(sum_Z2_U_err*sum_U_U_err-sum_Z_U_err*sum_UZ_U_err)/denom;
02018 double UZ=(sum_U_err*sum_UZ_U_err-sum_Z_U_err*sum_U_U_err)/denom;
02019
02022
02023 for(unsigned i=0; i<uu.size(); ++i){
02024 double uMean = U0+UZ*zz[i];
02025 chiUZind.push_back((pow((uMean-uu[i]),2))/e_Cxx[i]);
02026 chiUZ += (pow((uMean-uu[i]),2))/e_Cxx[i];
02027 }
02028 chiUZ = chiUZ/(uu.size()-2);
02029
02030 if(chiUZ>=chi2Norm_2D_){
02031 protoChiUCorrection = chiUZ/chi2Norm_2D_;
02032 for(unsigned i=0; i<uu.size(); ++i)
02033 e_Cxx[i]=e_Cxx[i]*protoChiUCorrection;
02034 }
02035
02037
02038 if(sqrt(protoChiUCorrection)>prePrunLimit_){
02039 chiContribution=max_element(chiUZind.begin(),chiUZind.end());
02040 maxContrIndex = chiContribution - chiUZind.begin();
02041
02042
02043
02044
02045
02046
02047
02048 }
02049
02050
02051 }
02052
02053 void CSCSegAlgoST::correctTheCovMatrix(CLHEP::HepMatrix &IC){
02054
02055 double condNumberCorr2=0.0;
02056 double detCov=0.0;
02057 double diag1=0.0;
02058 double diag2=0.0;
02059 double IC_12_corr=0.0;
02060 double IC_11_corr=0.0;
02061 if(!covToAnyNumberAll_){
02062
02063 condNumberCorr2=condSeed2_*IC(2,2);
02064 diag1=IC(1,1)*IC(2,2);
02065 diag2=IC(1,2)*IC(1,2);
02066 detCov=fabs(diag1-diag2);
02067 if((diag1<condNumberCorr2)&&(diag2<condNumberCorr2)){
02068 if(covToAnyNumber_)
02069 IC(1,2)=covAnyNumber_;
02070 else{
02071 IC_11_corr=condSeed1_+fabs(IC(1,2))/IC(2,2);
02072 IC(1,1)=IC_11_corr;
02073 }
02074 }
02075
02076 if(((detCov<condNumberCorr2)&&(diag1>condNumberCorr2))||
02077 ((diag2>condNumberCorr2)&&(detCov<condNumberCorr2)
02078 )){
02079 if(covToAnyNumber_)
02080 IC(1,2)=covAnyNumber_;
02081 else{
02082 IC_12_corr=sqrt(fabs(diag1-condNumberCorr2));
02083 if(IC(1,2)<0)
02084 IC(1,2)=(-1)*IC_12_corr;
02085 else
02086 IC(1,2)=IC_12_corr;
02087 }
02088 }
02089 }
02090 else{
02091 IC(1,2)=covAnyNumber_;
02092 }
02093 }
02094
02095 void CSCSegAlgoST::findDuplicates(std::vector<CSCSegment> & segments ){
02096
02097
02098
02099
02100
02101 for(std::vector<CSCSegment>::iterator it=segments.begin(); it != segments.end(); ++it) {
02102 std::vector<CSCSegment*> duplicateSegments;
02103 for(std::vector<CSCSegment>::iterator it2=segments.begin(); it2 != segments.end(); ++it2) {
02104
02105 bool allShared = true;
02106 if(it!=it2){
02107 allShared = it->sharesRecHits(*it2);
02108 }
02109 else{
02110 allShared = false;
02111 }
02112
02113 if(allShared){
02114 duplicateSegments.push_back(&(*it2));
02115 }
02116 }
02117 it->setDuplicateSegments(duplicateSegments);
02118 }
02119
02120 }
02121
02122