00001
00010 #include "CSCSegAlgoTC.h"
00011
00012 #include "DataFormats/CSCRecHit/interface/CSCSegment.h"
00013 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00014
00015
00016 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
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 CSCSegAlgoTC::CSCSegAlgoTC(const edm::ParameterSet& ps) : CSCSegmentAlgorithm(ps),
00030 myName("CSCSegAlgoTC") {
00031
00032 debugInfo = ps.getUntrackedParameter<bool>("verboseInfo");
00033
00034 dRPhiMax = ps.getParameter<double>("dRPhiMax");
00035 dPhiMax = ps.getParameter<double>("dPhiMax");
00036 dRPhiFineMax = ps.getParameter<double>("dRPhiFineMax");
00037 dPhiFineMax = ps.getParameter<double>("dPhiFineMax");
00038 chi2Max = ps.getParameter<double>("chi2Max");
00039 chi2ndfProbMin = ps.getParameter<double>("chi2ndfProbMin");
00040 minLayersApart = ps.getParameter<int>("minLayersApart");
00041 SegmentSorting = ps.getParameter<int>("SegmentSorting");
00042
00043 LogDebug("CSC") << myName << " has algorithm cuts set to: \n"
00044 << "--------------------------------------------------------------------\n"
00045 << "dRPhiMax = " << dRPhiMax << '\n'
00046 << "dPhiMax = " << dPhiMax << '\n'
00047 << "dRPhiFineMax = " << dRPhiFineMax << '\n'
00048 << "dPhiFineMax = " << dPhiFineMax << '\n'
00049 << "chi2Max = " << chi2Max << '\n'
00050 << "chi2ndfProbMin = " << chi2ndfProbMin << '\n'
00051 << "minLayersApart = " << minLayersApart << '\n'
00052 << "SegmentSorting = " << SegmentSorting << std::endl;
00053 }
00054
00055 std::vector<CSCSegment> CSCSegAlgoTC::run(const CSCChamber* aChamber, ChamberHitContainer rechits) {
00056 theChamber = aChamber;
00057 return buildSegments(rechits);
00058 }
00059
00060 std::vector<CSCSegment> CSCSegAlgoTC::buildSegments(ChamberHitContainer rechits) {
00061
00062
00063
00064 LogDebug("CSC") << "*********************************************";
00065 LogDebug("CSC") << "Start segment building in the new chamber: " << theChamber->specs()->chamberTypeName();
00066 LogDebug("CSC") << "*********************************************";
00067
00068 LayerIndex layerIndex(rechits.size());
00069
00070 for(unsigned int i = 0; i < rechits.size(); i++) {
00071
00072 layerIndex[i] = rechits[i]->cscDetId().layer();
00073 }
00074
00075 double z1 = theChamber->layer(1)->position().z();
00076 double z6 = theChamber->layer(6)->position().z();
00077
00078 if ( z1 > 0. ) {
00079 if ( z1 > z6 ) {
00080 reverse(layerIndex.begin(), layerIndex.end());
00081 reverse(rechits.begin(), rechits.end());
00082 }
00083 }
00084 else if ( z1 < 0. ) {
00085 if ( z1 < z6 ) {
00086 reverse(layerIndex.begin(), layerIndex.end());
00087 reverse(rechits.begin(), rechits.end());
00088 }
00089 }
00090
00091 if (debugInfo) {
00092
00093 dumpHits(rechits);
00094 }
00095
00096 if (rechits.size() < 2) {
00097 LogDebug("CSC") << myName << ": " << rechits.size() <<
00098 " hit(s) in chamber is not enough to build a segment.\n";
00099 return std::vector<CSCSegment>();
00100 }
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119 std::vector<CSCSegment> segments;
00120
00121 ChamberHitContainerCIt ib = rechits.begin();
00122 ChamberHitContainerCIt ie = rechits.end();
00123
00124 for (ChamberHitContainerCIt i1 = ib; i1 != ie; ++i1) {
00125
00126 int layer1 = layerIndex[i1-ib];
00127 const CSCRecHit2D* h1 = *i1;
00128
00129 for (ChamberHitContainerCIt i2 = ie-1; i2 != i1; --i2) {
00130
00131 int layer2 = layerIndex[i2-ib];
00132
00133 if (abs(layer2 - layer1) < minLayersApart)
00134 break;
00135
00136 const CSCRecHit2D* h2 = *i2;
00137
00138 if (areHitsCloseInLocalX(h1, h2) && areHitsCloseInGlobalPhi(h1, h2)) {
00139
00140 proto_segment.clear();
00141
00142 const CSCLayer* l1 = theChamber->layer(h1->cscDetId().layer());
00143 GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
00144 const CSCLayer* l2 = theChamber->layer(h2->cscDetId().layer());
00145 GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
00146 LogDebug("CSC") << "start new segment from hits " << "h1: " << gp1 << " - h2: " << gp2 << "\n";
00147
00148 if (!addHit(h1, layer1)) {
00149 LogDebug("CSC") << " failed to add hit h1\n";
00150 continue;
00151 }
00152
00153 if (!addHit(h2, layer2)) {
00154 LogDebug("CSC") << " failed to add hit h2\n";
00155 continue;
00156 }
00157
00158 tryAddingHitsToSegment(rechits, i1, i2);
00159
00160
00161 if (proto_segment.empty()) {
00162
00163 LogDebug("CSC") << "No segment has been found !!!\n";
00164 }
00165 else {
00166
00167
00168 AlgebraicSymMatrix error_matrix = calculateError();
00169
00170
00171
00172 flipErrors( error_matrix );
00173
00174 candidates.push_back(proto_segment);
00175 origins.push_back(theOrigin);
00176 directions.push_back(theDirection);
00177 errors.push_back(error_matrix);
00178 chi2s.push_back(theChi2);
00179 LogDebug("CSC") << "Found a segment !!!\n";
00180 }
00181 }
00182 }
00183 }
00184
00185
00186 pruneTheSegments(rechits);
00187
00188
00189 for(unsigned int i=0; i < candidates.size(); i++) {
00190
00191 CSCSegment temp(candidates[i], origins[i], directions[i], errors[i], chi2s[i]);
00192 segments.push_back(temp);
00193 }
00194
00195 candidates.clear();
00196 origins.clear();
00197 directions.clear();
00198 errors.clear();
00199 chi2s.clear();
00200
00201
00202 return segments;
00203 }
00204
00205 void CSCSegAlgoTC::tryAddingHitsToSegment(const ChamberHitContainer& rechits,
00206 const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2) {
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218 ChamberHitContainerCIt ib = rechits.begin();
00219 ChamberHitContainerCIt ie = rechits.end();
00220
00221 for (ChamberHitContainerCIt i = ib; i != ie; ++i) {
00222
00223 if ( i == i1 || i == i2 )
00224 continue;
00225
00226 int layer = (*i)->cscDetId().layer();
00227 const CSCRecHit2D* h = *i;
00228
00229 if (isHitNearSegment(h)) {
00230
00231 const CSCLayer* l1 = theChamber->layer(h->cscDetId().layer());
00232 GlobalPoint gp1 = l1->toGlobal(h->localPosition());
00233 LogDebug("CSC") << " hit at global " << gp1 << " is near segment\n.";
00234
00235 if (hasHitOnLayer(layer)) {
00236 if (proto_segment.size() <= 2) {
00237 LogDebug("CSC") << " " << proto_segment.size()
00238 << " hits on segment...skip hit on same layer.\n";
00239 continue;
00240 }
00241
00242 compareProtoSegment(h, layer);
00243 }
00244 else
00245 increaseProtoSegment(h, layer);
00246 }
00247 }
00248 }
00249
00250 bool CSCSegAlgoTC::addHit(const CSCRecHit2D* aHit, int layer) {
00251
00252
00253
00254
00255 bool ok = true;
00256 ChamberHitContainer::const_iterator it;
00257
00258 for(it = proto_segment.begin(); it != proto_segment.end(); it++)
00259 if (((*it)->cscDetId().layer() == layer) && (aHit != *it))
00260 return false;
00261
00262 if (ok) {
00263 proto_segment.push_back(aHit);
00264 updateParameters();
00265 }
00266 return ok;
00267 }
00268
00269 bool CSCSegAlgoTC::replaceHit(const CSCRecHit2D* h, int layer) {
00270
00271
00272 ChamberHitContainer::iterator it;
00273 for (it = proto_segment.begin(); it != proto_segment.end();) {
00274 if ((*it)->cscDetId().layer() == layer)
00275 it = proto_segment.erase(it);
00276 else
00277 ++it;
00278 }
00279
00280 return addHit(h, layer);
00281 }
00282
00283 void CSCSegAlgoTC::updateParameters() {
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295 int nh = proto_segment.size();
00296
00297
00298 if (nh < 2)
00299 return;
00300
00301 if (nh == 2) {
00302
00303
00304
00305 ChamberHitContainer::const_iterator ih = proto_segment.begin();
00306 int il1 = (*ih)->cscDetId().layer();
00307 const CSCRecHit2D& h1 = (**ih);
00308 ++ih;
00309 int il2 = (*ih)->cscDetId().layer();
00310 const CSCRecHit2D& h2 = (**ih);
00311
00312
00313 if (il1 == il2)
00314 return;
00315
00316 const CSCLayer* layer1 = theChamber->layer(il1);
00317 const CSCLayer* layer2 = theChamber->layer(il2);
00318
00319 GlobalPoint h1glopos = layer1->toGlobal(h1.localPosition());
00320 GlobalPoint h2glopos = layer2->toGlobal(h2.localPosition());
00321
00322
00323 theOrigin = h1.localPosition();
00324
00325
00326 LocalPoint h1pos = theChamber->toLocal(h1glopos);
00327 LocalPoint h2pos = theChamber->toLocal(h2glopos);
00328
00329 float dz = h2pos.z()-h1pos.z();
00330 uz = (h2pos.x()-h1pos.x())/dz ;
00331 vz = (h2pos.y()-h1pos.y())/dz ;
00332
00333 theChi2 = 0.;
00334 }
00335 else if (nh > 2) {
00336
00337
00338 fitSlopes();
00339 fillChiSquared();
00340 }
00341
00342 fillLocalDirection();
00343 }
00344
00345 void CSCSegAlgoTC::fitSlopes() {
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426 CLHEP::HepMatrix M(4,4,0);
00427 CLHEP::HepVector B(4,0);
00428
00429 ChamberHitContainer::const_iterator ih = proto_segment.begin();
00430
00431 for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
00432
00433 const CSCRecHit2D& hit = (**ih);
00434 const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00435 GlobalPoint gp = layer->toGlobal(hit.localPosition());
00436 LocalPoint lp = theChamber->toLocal(gp);
00437
00438
00439 double u = lp.x();
00440 double v = lp.y();
00441 double z = lp.z();
00442
00443
00444 CLHEP::HepMatrix IC(2,2);
00445 IC(1,1) = hit.localPositionError().xx();
00446 IC(1,2) = hit.localPositionError().xy();
00447 IC(2,1) = IC(1,2);
00448 IC(2,2) = hit.localPositionError().yy();
00449
00450
00451 int ierr = 0;
00452 IC.invert(ierr);
00453 if (ierr != 0) {
00454 LogDebug("CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n";
00455
00456
00457 }
00458
00459
00460 M(1,1) += IC(1,1);
00461 M(1,2) += IC(1,2);
00462 M(1,3) += IC(1,1) * z;
00463 M(1,4) += IC(1,2) * z;
00464 B(1) += u * IC(1,1) + v * IC(1,2);
00465
00466 M(2,1) += IC(2,1);
00467 M(2,2) += IC(2,2);
00468 M(2,3) += IC(2,1) * z;
00469 M(2,4) += IC(2,2) * z;
00470 B(2) += u * IC(2,1) + v * IC(2,2);
00471
00472 M(3,1) += IC(1,1) * z;
00473 M(3,2) += IC(1,2) * z;
00474 M(3,3) += IC(1,1) * z * z;
00475 M(3,4) += IC(1,2) * z * z;
00476 B(3) += ( u * IC(1,1) + v * IC(1,2) ) * z;
00477
00478 M(4,1) += IC(2,1) * z;
00479 M(4,2) += IC(2,2) * z;
00480 M(4,3) += IC(2,1) * z * z;
00481 M(4,4) += IC(2,2) * z * z;
00482 B(4) += ( u * IC(2,1) + v * IC(2,2) ) * z;
00483 }
00484
00485
00486
00487 CLHEP::HepVector p = solve(M, B);
00488
00489
00490 theOrigin = LocalPoint(p(1), p(2), 0.);
00491 uz = p(3);
00492 vz = p(4);
00493 }
00494
00495 void CSCSegAlgoTC::fillChiSquared() {
00496
00497
00498
00499
00500
00501 float u0 = theOrigin.x();
00502 float v0 = theOrigin.y();
00503 double chsq = 0.;
00504
00505 ChamberHitContainer::const_iterator ih;
00506 for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
00507
00508 const CSCRecHit2D& hit = (**ih);
00509 const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00510 GlobalPoint gp = layer->toGlobal(hit.localPosition());
00511 LocalPoint lp = theChamber->toLocal(gp);
00512
00513 double hu = lp.x();
00514 double hv = lp.y();
00515 double hz = lp.z();
00516
00517 double du = u0 + uz * hz - hu;
00518 double dv = v0 + vz * hz - hv;
00519
00520 CLHEP::HepMatrix IC(2,2);
00521 IC(1,1) = hit.localPositionError().xx();
00522 IC(1,2) = hit.localPositionError().xy();
00523 IC(2,1) = IC(1,2);
00524 IC(2,2) = hit.localPositionError().yy();
00525
00526
00527 int ierr = 0;
00528 IC.invert(ierr);
00529 if (ierr != 0) {
00530 LogDebug("CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";
00531
00532
00533 }
00534
00535 chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
00536 }
00537 theChi2 = chsq;
00538 }
00539
00540 void CSCSegAlgoTC::fillLocalDirection() {
00541
00542
00543
00544
00545 double dxdz = uz;
00546 double dydz = vz;
00547 double dz = 1./sqrt(1. + dxdz*dxdz + dydz*dydz);
00548 double dx = dz*dxdz;
00549 double dy = dz*dydz;
00550 LocalVector localDir(dx,dy,dz);
00551
00552
00553
00554
00555
00556 double globalZpos = ( theChamber->toGlobal( theOrigin ) ).z();
00557 double globalZdir = ( theChamber->toGlobal( localDir ) ).z();
00558 double directionSign = globalZpos * globalZdir;
00559
00560 theDirection = (directionSign * localDir).unit();
00561 }
00562
00563 float CSCSegAlgoTC::phiAtZ(float z) const {
00564
00565
00566 const CSCLayer* l1 = theChamber->layer(1);
00567 GlobalPoint gp = l1->toGlobal(theOrigin);
00568 GlobalVector gv = l1->toGlobal(theDirection);
00569
00570 float x = gp.x() + (gv.x()/gv.z())*(z - gp.z());
00571 float y = gp.y() + (gv.y()/gv.z())*(z - gp.z());
00572 float phi = atan2(y, x);
00573 if (phi < 0.f)
00574 phi += 2. * M_PI;
00575
00576 return phi ;
00577 }
00578
00579 void CSCSegAlgoTC::compareProtoSegment(const CSCRecHit2D* h, int layer) {
00580
00581
00582 double oldChi2 = theChi2;
00583 LocalPoint oldOrigin = theOrigin;
00584 LocalVector oldDirection = theDirection;
00585 ChamberHitContainer oldSegment = proto_segment;
00586
00587 bool ok = replaceHit(h, layer);
00588
00589 if (ok) {
00590 LogDebug("CSC") << " hit in same layer as a hit on segment; try replacing old one..."
00591 << " chi2 new: " << theChi2 << " old: " << oldChi2 << "\n";
00592 }
00593
00594 if ((theChi2 < oldChi2) && (ok)) {
00595 LogDebug("CSC") << " segment with replaced hit is better.\n";
00596 }
00597 else {
00598 proto_segment = oldSegment;
00599 theChi2 = oldChi2;
00600 theOrigin = oldOrigin;
00601 theDirection = oldDirection;
00602 }
00603 }
00604
00605 void CSCSegAlgoTC::increaseProtoSegment(const CSCRecHit2D* h, int layer) {
00606
00607 double oldChi2 = theChi2;
00608 LocalPoint oldOrigin = theOrigin;
00609 LocalVector oldDirection = theDirection;
00610 ChamberHitContainer oldSegment = proto_segment;
00611
00612 bool ok = addHit(h, layer);
00613
00614 if (ok) {
00615 LogDebug("CSC") << " hit in new layer: added to segment, new chi2: "
00616 << theChi2 << "\n";
00617 }
00618
00619 int ndf = 2*proto_segment.size() - 4;
00620
00621 if (ok && ((ndf <= 0) || (theChi2/ndf < chi2Max))) {
00622 LogDebug("CSC") << " segment with added hit is good.\n" ;
00623 }
00624 else {
00625 proto_segment = oldSegment;
00626 theChi2 = oldChi2;
00627 theOrigin = oldOrigin;
00628 theDirection = oldDirection;
00629 }
00630 }
00631
00632 bool CSCSegAlgoTC::areHitsCloseInLocalX(const CSCRecHit2D* h1, const CSCRecHit2D* h2) const {
00633 float h1x = h1->localPosition().x();
00634 float h2x = h2->localPosition().x();
00635 float deltaX = (h1->localPosition()-h2->localPosition()).x();
00636 LogDebug("CSC") << " Hits at local x= " << h1x << ", "
00637 << h2x << " have separation= " << deltaX;
00638 return (fabs(deltaX) < (dRPhiMax))? true:false;
00639 }
00640
00641 bool CSCSegAlgoTC::areHitsCloseInGlobalPhi(const CSCRecHit2D* h1, const CSCRecHit2D* h2) const {
00642
00643 const CSCLayer* l1 = theChamber->layer(h1->cscDetId().layer());
00644 GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
00645 const CSCLayer* l2 = theChamber->layer(h2->cscDetId().layer());
00646 GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
00647
00648 float h1p = gp1.phi();
00649 float h2p = gp2.phi();
00650 float dphi12 = h1p - h2p;
00651
00652
00653 if (dphi12 < -M_PI)
00654 dphi12 += 2.*M_PI;
00655 if (dphi12 > M_PI)
00656 dphi12 -= 2.*M_PI;
00657 LogDebug("CSC") << " Hits at global phi= " << h1p << ", "
00658 << h2p << " have separation= " << dphi12;
00659 return (fabs(dphi12) < dPhiMax)? true:false;
00660 }
00661
00662 bool CSCSegAlgoTC::isHitNearSegment(const CSCRecHit2D* h) const {
00663
00664
00665
00666
00667
00668
00669
00670 const CSCLayer* l1 = theChamber->layer(h->cscDetId().layer());
00671 GlobalPoint hp = l1->toGlobal(h->localPosition());
00672
00673 float hphi = hp.phi();
00674 if (hphi < 0.)
00675 hphi += 2.*M_PI;
00676 float sphi = phiAtZ(hp.z());
00677 float phidif = sphi-hphi;
00678 if (phidif < 0.)
00679 phidif += 2.*M_PI;
00680 if (phidif > M_PI)
00681 phidif -= 2.*M_PI;
00682
00683 float dRPhi = fabs(phidif)*hp.perp();
00684 LogDebug("CSC") << " is hit at phi_h= " << hphi << " near segment phi_seg= " << sphi
00685 << "? is " << dRPhi << "<" << dRPhiFineMax << " ? "
00686 << " and is |" << phidif << "|<" << dPhiFineMax << " ?";
00687
00688 return ((dRPhi < dRPhiFineMax) &&
00689 (fabs(phidif) < dPhiFineMax))? true:false;
00690 }
00691
00692 bool CSCSegAlgoTC::hasHitOnLayer(int layer) const {
00693
00694
00695 ChamberHitContainerCIt it;
00696
00697 for(it = proto_segment.begin(); it != proto_segment.end(); it++)
00698 if ((*it)->cscDetId().layer() == layer)
00699 return true;
00700
00701 return false;
00702 }
00703
00704 void CSCSegAlgoTC::dumpHits(const ChamberHitContainer& rechits) const {
00705
00706
00707 ChamberHitContainerCIt it;
00708
00709 for(it=rechits.begin(); it!=rechits.end(); it++) {
00710
00711 const CSCLayer* l1 = theChamber->layer((*it)->cscDetId().layer());
00712 GlobalPoint gp1 = l1->toGlobal((*it)->localPosition());
00713
00714 LogDebug("CSC") << "Global pos.: " << gp1 << ", phi: " << gp1.phi() << ". Local position: "
00715 << (*it)->localPosition() << ", phi: "
00716 << (*it)->localPosition().phi() << ". Layer: "
00717 << (*it)->cscDetId().layer() << "\n";
00718 }
00719 }
00720
00721
00722 bool CSCSegAlgoTC::isSegmentGood(std::vector<ChamberHitContainer>::iterator seg, std::vector<double>::iterator chi2,
00723 const ChamberHitContainer& rechitsInChamber, BoolContainer& used) const {
00724
00725
00726
00727
00728
00729
00730
00731
00732
00733 unsigned int iadd = (rechitsInChamber.size() > 20 )? 1 : 0;
00734
00735 if (seg->size() < 3 + iadd)
00736 return false;
00737
00738
00739
00740
00741
00742
00743 if( SegmentSorting == 2 ){
00744 if( (*chi2) != 0 && ((2*seg->size())-4) >0 ) {
00745 if ( ChiSquaredProbability((*chi2),(double)(2*seg->size()-4)) < chi2ndfProbMin ) {
00746 return false;
00747 }
00748 }
00749 if((*chi2) == 0 ) return false;
00750 }
00751
00752
00753 for(unsigned int ish = 0; ish < seg->size(); ++ish) {
00754
00755 ChamberHitContainerCIt ib = rechitsInChamber.begin();
00756
00757 for(ChamberHitContainerCIt ic = ib; ic != rechitsInChamber.end(); ++ic) {
00758
00759 if(((*seg)[ish] == (*ic)) && used[ic-ib])
00760 return false;
00761 }
00762 }
00763
00764 return true;
00765 }
00766
00767 void CSCSegAlgoTC::flagHitsAsUsed(std::vector<ChamberHitContainer>::iterator seg, const ChamberHitContainer& rechitsInChamber,
00768 BoolContainer& used) const {
00769
00770
00771
00772 ChamberHitContainerCIt ib = rechitsInChamber.begin();
00773
00774 for(unsigned int ish = 0; ish < seg->size(); ish++) {
00775
00776 for(ChamberHitContainerCIt iu = ib; iu != rechitsInChamber.end(); ++iu)
00777 if((*seg)[ish] == (*iu))
00778 used[iu-ib] = true;
00779 }
00780 }
00781
00782 void CSCSegAlgoTC::pruneTheSegments(const ChamberHitContainer& rechitsInChamber) {
00783
00784
00785
00786
00787 if (candidates.empty())
00788 return;
00789
00790
00791 BoolContainer used(rechitsInChamber.size(), false);
00792
00793
00794 segmentSort();
00795
00796
00797
00798
00799
00800 std::vector<ChamberHitContainer>::iterator is;
00801 std::vector<double>::iterator ichi = chi2s.begin();
00802 std::vector<AlgebraicSymMatrix>::iterator iErr = errors.begin();
00803 std::vector<LocalPoint>::iterator iOrig = origins.begin();
00804 std::vector<LocalVector>::iterator iDir = directions.begin();
00805
00806 for (is = candidates.begin(); is != candidates.end(); ) {
00807
00808 bool goodSegment = isSegmentGood(is, ichi, rechitsInChamber, used);
00809
00810 if (goodSegment) {
00811 LogDebug("CSC") << "Accepting segment: ";
00812
00813 flagHitsAsUsed(is, rechitsInChamber, used);
00814 ++is;
00815 ++ichi;
00816 ++iErr;
00817 ++iOrig;
00818 ++iDir;
00819 }
00820 else {
00821 LogDebug("CSC") << "Rejecting segment: ";
00822 is = candidates.erase(is);
00823 ichi = chi2s.erase(ichi);
00824 iErr = errors.erase(iErr);
00825 iOrig = origins.erase(iOrig);
00826 iDir = directions.erase(iDir);
00827 }
00828 }
00829 }
00830
00831 void CSCSegAlgoTC::segmentSort() {
00832
00833
00834
00835 for(unsigned int i=0; i<candidates.size()-1; i++) {
00836 for(unsigned int j=i+1; j<candidates.size(); j++) {
00837
00838 ChamberHitContainer s1 = candidates[i];
00839 ChamberHitContainer s2 = candidates[j];
00840 if (i == j)
00841 continue;
00842
00843 int n1 = candidates[i].size();
00844 int n2 = candidates[j].size();
00845
00846 if( SegmentSorting == 2 ){
00847 if ( n2 > n1 ) {
00848 ChamberHitContainer temp = candidates[j];
00849 candidates[j] = candidates[i];
00850 candidates[i] = temp;
00851
00852 double temp1 = chi2s[j];
00853 chi2s[j] = chi2s[i];
00854 chi2s[i] = temp1;
00855
00856 AlgebraicSymMatrix temp2 = errors[j];
00857 errors[j] = errors[i];
00858 errors[i] = temp2;
00859
00860 LocalPoint temp3 = origins[j];
00861 origins[j] = origins[i];
00862 origins[i] = temp3;
00863
00864 LocalVector temp4 = directions[j];
00865 directions[j] = directions[i];
00866 directions[i] = temp4;
00867 }
00868
00869 if(chi2s[i] != 0. && 2*n2-4 > 0 ) {
00870 if( n2 == n1 && (ChiSquaredProbability( chi2s[i],(double)(2*n1-4)) < ChiSquaredProbability(chi2s[j],(double)(2*n2-4))) ){
00871 ChamberHitContainer temp = candidates[j];
00872 candidates[j] = candidates[i];
00873 candidates[i] = temp;
00874
00875 double temp1 = chi2s[j];
00876 chi2s[j] = chi2s[i];
00877 chi2s[i] = temp1;
00878
00879 AlgebraicSymMatrix temp2 = errors[j];
00880 errors[j] = errors[i];
00881 errors[i] = temp2;
00882
00883 LocalPoint temp3 = origins[j];
00884 origins[j] = origins[i];
00885 origins[i] = temp3;
00886
00887 LocalVector temp4 = directions[j];
00888 directions[j] = directions[i];
00889 directions[i] = temp4;
00890 }
00891 }
00892 }
00893 else if( SegmentSorting == 1 ){
00894 if ((chi2s[i]/n1) > (chi2s[j]/n2)) {
00895
00896 ChamberHitContainer temp = candidates[j];
00897 candidates[j] = candidates[i];
00898 candidates[i] = temp;
00899
00900 double temp1 = chi2s[j];
00901 chi2s[j] = chi2s[i];
00902 chi2s[i] = temp1;
00903
00904 AlgebraicSymMatrix temp2 = errors[j];
00905 errors[j] = errors[i];
00906 errors[i] = temp2;
00907
00908 LocalPoint temp3 = origins[j];
00909 origins[j] = origins[i];
00910 origins[i] = temp3;
00911
00912 LocalVector temp4 = directions[j];
00913 directions[j] = directions[i];
00914 directions[i] = temp4;
00915 }
00916 }
00917 else{
00918 LogDebug("CSC") << "No valid segment sorting specified - BAD !!!\n";
00919 }
00920 }
00921 }
00922 }
00923
00924 AlgebraicSymMatrix CSCSegAlgoTC::calculateError() const {
00925
00926 AlgebraicSymMatrix weights = weightMatrix();
00927 AlgebraicMatrix A = derivativeMatrix();
00928
00929
00930
00931 int ierr;
00932 AlgebraicSymMatrix result = weights.similarityT(A);
00933 result.invert(ierr);
00934
00935
00936 return result;
00937 }
00938
00939 CLHEP::HepMatrix CSCSegAlgoTC::derivativeMatrix() const {
00940
00941 ChamberHitContainer::const_iterator it;
00942 int nhits = proto_segment.size();
00943 CLHEP::HepMatrix matrix(2*nhits, 4);
00944 int row = 0;
00945
00946 for(it = proto_segment.begin(); it != proto_segment.end(); ++it) {
00947
00948 const CSCRecHit2D& hit = (**it);
00949 const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00950 GlobalPoint gp = layer->toGlobal(hit.localPosition());
00951 LocalPoint lp = theChamber->toLocal(gp);
00952 float z = lp.z();
00953 ++row;
00954 matrix(row, 1) = 1.;
00955 matrix(row, 3) = z;
00956 ++row;
00957 matrix(row, 2) = 1.;
00958 matrix(row, 4) = z;
00959 }
00960 return matrix;
00961 }
00962
00963
00964 AlgebraicSymMatrix CSCSegAlgoTC::weightMatrix() const {
00965
00966 std::vector<const CSCRecHit2D*>::const_iterator it;
00967 int nhits = proto_segment.size();
00968 AlgebraicSymMatrix matrix(2*nhits, 0);
00969 int row = 0;
00970
00971 for (it = proto_segment.begin(); it != proto_segment.end(); ++it) {
00972
00973 const CSCRecHit2D& hit = (**it);
00974 ++row;
00975 matrix(row, row) = hit.localPositionError().xx();
00976 matrix(row, row+1) = hit.localPositionError().xy();
00977 ++row;
00978 matrix(row, row-1) = hit.localPositionError().xy();
00979 matrix(row, row) = hit.localPositionError().yy();
00980 }
00981 int ierr;
00982 matrix.invert(ierr);
00983 return matrix;
00984 }
00985
00986 void CSCSegAlgoTC::flipErrors( AlgebraicSymMatrix& a ) const {
00987
00988
00989
00990 AlgebraicSymMatrix hold( a );
00991
00992
00993 a(1,1) = hold(3,3);
00994 a(1,2) = hold(3,4);
00995 a(2,1) = hold(4,3);
00996 a(2,2) = hold(4,4);
00997
00998
00999 a(3,3) = hold(1,1);
01000 a(3,4) = hold(1,2);
01001 a(4,3) = hold(2,1);
01002 a(4,4) = hold(2,2);
01003
01004
01005
01006 }
01007