Go to the documentation of this file.00001
00007 #include <RecoMuon/MuonSeedGenerator/src/MuonSeedCleaner.h>
00008
00009
00010 #include <DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h>
00011 #include <DataFormats/TrajectorySeed/interface/TrajectorySeed.h>
00012 #include <DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h>
00013 #include <DataFormats/CSCRecHit/interface/CSCRecHit2D.h>
00014 #include <DataFormats/CSCRecHit/interface/CSCSegmentCollection.h>
00015 #include <DataFormats/CSCRecHit/interface/CSCSegment.h>
00016 #include <DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h>
00017 #include <DataFormats/DTRecHit/interface/DTRecSegment4D.h>
00018
00019
00020 #include <Geometry/CommonDetUnit/interface/GeomDet.h>
00021 #include <TrackingTools/DetLayers/interface/DetLayer.h>
00022 #include <RecoMuon/MeasurementDet/interface/MuonDetLayerMeasurements.h>
00023 #include <RecoMuon/DetLayers/interface/MuonDetLayerGeometry.h>
00024 #include <RecoMuon/Records/interface/MuonRecoGeometryRecord.h>
00025
00026
00027 #include <TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h>
00028 #include <DataFormats/TrajectoryState/interface/PTrajectoryStateOnDet.h>
00029
00030
00031 #include <FWCore/ParameterSet/interface/ParameterSet.h>
00032 #include "FWCore/Framework/interface/Event.h"
00033 #include <FWCore/Framework/interface/ESHandle.h>
00034 #include <FWCore/MessageLogger/interface/MessageLogger.h>
00035 #include <DataFormats/Common/interface/Handle.h>
00036
00037
00038 #include <vector>
00039 #include <deque>
00040 #include <utility>
00041
00042
00043
00044 static bool lengthSorting(const TrajectorySeed s1, const TrajectorySeed s2) { return ( s1.nHits() > s2.nHits() ); }
00045
00046
00047
00048
00049 MuonSeedCleaner::MuonSeedCleaner(const edm::ParameterSet& pset){
00050
00051
00052 debug = pset.getParameter<bool>("DebugMuonSeed");
00053
00054
00055 edm::ParameterSet serviceParameters = pset.getParameter<edm::ParameterSet>("ServiceParameters");
00056 theService = new MuonServiceProxy(serviceParameters);
00057
00058 }
00059
00060
00061
00062
00063 MuonSeedCleaner::~MuonSeedCleaner(){
00064
00065 if (theService) delete theService;
00066 }
00067
00068
00069
00070
00071
00072
00073
00074 std::vector<TrajectorySeed> MuonSeedCleaner::seedCleaner(const edm::EventSetup& eventSetup, std::vector<TrajectorySeed>& seeds ) {
00075
00076 theService->update(eventSetup);
00077
00078 std::vector<TrajectorySeed> FinalSeeds;
00079
00080
00081 std::vector<SeedContainer> theCollection = GroupSeeds(seeds);
00082
00083
00084 for (size_t i=0; i< theCollection.size(); i++ ) {
00085
00086
00087 SeedContainer goodSeeds = SeedCandidates( theCollection[i], true );
00088 SeedContainer otherSeeds = SeedCandidates( theCollection[i], false );
00089 if ( MomentumFilter( goodSeeds ) ) {
00090
00091 TrajectorySeed bestSeed = Chi2LengthSelection( goodSeeds );
00092 FinalSeeds.push_back( bestSeed );
00093
00094 GlobalPoint seedgp = SeedPosition( bestSeed );
00095 double eta = fabs( seedgp.eta() );
00096 if ( goodSeeds.size() > 2 && eta > 1.5 ) {
00097 TrajectorySeed anotherSeed = MoreRecHits( goodSeeds );
00098 FinalSeeds.push_back( anotherSeed );
00099 }
00100 }
00101 else if( MomentumFilter( otherSeeds ) ) {
00102
00103 TrajectorySeed bestSeed = MoreRecHits( otherSeeds );
00104 FinalSeeds.push_back( bestSeed );
00105
00106 GlobalPoint seedgp = SeedPosition( bestSeed );
00107 double eta = fabs( seedgp.eta() );
00108 if ( otherSeeds.size() > 2 && eta > 1.5 ) {
00109 TrajectorySeed anotherSeed = LeanHighMomentum( otherSeeds );
00110 FinalSeeds.push_back( anotherSeed );
00111 }
00112 }
00113 else {
00114
00115 TrajectorySeed bestSeed = LeanHighMomentum( theCollection[i] );
00116 FinalSeeds.push_back( bestSeed );
00117
00118 GlobalPoint seedgp = SeedPosition( bestSeed );
00119 double eta = fabs( seedgp.eta() );
00120 if ( theCollection.size() > 2 && eta > 1.5 ) {
00121 TrajectorySeed anotherSeed = BiggerCone( theCollection[i] );
00122 FinalSeeds.push_back( anotherSeed );
00123 }
00124 }
00125 }
00126 return FinalSeeds ;
00127
00128 }
00129
00130
00131 TrajectorySeed MuonSeedCleaner::Chi2LengthSelection(std::vector<TrajectorySeed>& seeds ) {
00132
00133 if ( seeds.size() == 1 ) return seeds[0];
00134
00135 int winner = 0 ;
00136 int moreHits = 0 ;
00137 double bestChi2 = 99999.;
00138 for ( size_t i = 0; i < seeds.size(); i++ ) {
00139
00140
00141
00142
00143
00144
00145 double theChi2 = SeedChi2( seeds[i] );
00146 double dChi2 = fabs( 1. - (theChi2 / bestChi2)) ;
00147 int theHits = seeds[i].nHits() ;
00148 int dHits = theHits - moreHits ;
00149
00150
00151
00152 if ( theChi2 < bestChi2 && dChi2 > 0.05 ) {
00153 winner = static_cast<int>(i) ;
00154 bestChi2 = theChi2 ;
00155 moreHits = theHits ;
00156 }
00157
00158 if ( theChi2 >= bestChi2 && dChi2 < 0.05 && dHits > 0 ) {
00159 winner = static_cast<int>(i) ;
00160 bestChi2 = theChi2 ;
00161 moreHits = theHits ;
00162 }
00163
00164 }
00165
00166 TrajectorySeed theSeed = seeds[winner];
00167 seeds.erase( seeds.begin()+winner );
00168 return theSeed;
00169 }
00170
00171 TrajectorySeed MuonSeedCleaner::BiggerCone(std::vector<TrajectorySeed>& seeds ) {
00172
00173 if ( seeds.size() == 1 ) return seeds[0];
00174
00175 float biggerProjErr = 9999.;
00176 int winner = 0 ;
00177 AlgebraicSymMatrix mat(5,0) ;
00178 for ( size_t i = 0; i < seeds.size(); i++ ) {
00179
00180 edm::OwnVector<TrackingRecHit>::const_iterator r1 = seeds[i].recHits().first ;
00181 mat = r1->parametersError().similarityT( r1->projectionMatrix() );
00182
00183 int NRecHits = NRecHitsFromSegment( *r1 );
00184
00185 float ddx = mat[1][1];
00186 float ddy = mat[2][2];
00187 float dxx = mat[3][3];
00188 float dyy = mat[4][4];
00189 float projectErr = sqrt( (ddx*10000.) + (ddy*10000.) + dxx + dyy ) ;
00190
00191 if ( NRecHits < 5 ) continue ;
00192 if ( projectErr < biggerProjErr ) continue;
00193
00194 winner = static_cast<int>(i) ;
00195 biggerProjErr = projectErr ;
00196 }
00197 TrajectorySeed theSeed = seeds[winner];
00198 seeds.erase( seeds.begin()+winner );
00199 return theSeed;
00200
00201 }
00202
00203 TrajectorySeed MuonSeedCleaner::LeanHighMomentum( std::vector<TrajectorySeed>& seeds ) {
00204
00205 if ( seeds.size() == 1 ) return seeds[0];
00206
00207 double highestPt = 0. ;
00208 int winner = 0 ;
00209 for ( size_t i = 0; i < seeds.size(); i++ ) {
00210 GlobalVector mom = SeedMomentum( seeds[i] );
00211 double pt = sqrt( (mom.x()*mom.x()) + (mom.y()*mom.y()) );
00212 if ( pt > highestPt ) {
00213 winner = static_cast<int>(i) ;
00214 highestPt = pt ;
00215 }
00216 }
00217 TrajectorySeed theSeed = seeds[winner];
00218 seeds.erase( seeds.begin()+winner );
00219 return theSeed;
00220
00221 }
00222
00223
00224 TrajectorySeed MuonSeedCleaner::MoreRecHits(std::vector<TrajectorySeed>& seeds ) {
00225
00226 if ( seeds.size() == 1 ) return seeds[0];
00227
00228 int winner = 0 ;
00229 int moreHits = 0 ;
00230 double betterChi2 = 99999.;
00231 for ( size_t i = 0; i < seeds.size(); i++ ) {
00232
00233 int theHits = 0;
00234 for (edm::OwnVector<TrackingRecHit>::const_iterator r1 = seeds[i].recHits().first; r1 != seeds[i].recHits().second; r1++){
00235 theHits += NRecHitsFromSegment( *r1 );
00236 }
00237
00238 double theChi2 = SeedChi2( seeds[i] );
00239
00240 if ( theHits == moreHits && theChi2 < betterChi2 ) {
00241 betterChi2 = theChi2;
00242 winner = static_cast<int>(i) ;
00243 }
00244 if ( theHits > moreHits ) {
00245 moreHits = theHits ;
00246 betterChi2 = theChi2;
00247 winner = static_cast<int>(i) ;
00248 }
00249 }
00250 TrajectorySeed theSeed = seeds[winner];
00251 seeds.erase( seeds.begin()+winner );
00252 return theSeed;
00253 }
00254
00255
00256 SeedContainer MuonSeedCleaner::LengthFilter(std::vector<TrajectorySeed>& seeds ) {
00257
00258 SeedContainer longSeeds;
00259 int NSegs = 0;
00260 for (size_t i = 0; i< seeds.size(); i++) {
00261
00262 int theLength = static_cast<int>( seeds[i].nHits());
00263 if ( theLength > NSegs ) {
00264 NSegs = theLength ;
00265 longSeeds.clear();
00266 longSeeds.push_back( seeds[i] );
00267 }
00268 else if ( theLength == NSegs ) {
00269 longSeeds.push_back( seeds[i] );
00270 } else {
00271 continue;
00272 }
00273 }
00274
00275
00276 return longSeeds ;
00277
00278 }
00279
00280
00281 bool MuonSeedCleaner::MomentumFilter(std::vector<TrajectorySeed>& seeds ) {
00282
00283 bool findgoodMomentum = false;
00284 SeedContainer goodMomentumSeeds = seeds;
00285 seeds.clear();
00286 for ( size_t i = 0; i < goodMomentumSeeds.size(); i++ ) {
00287 GlobalVector mom = SeedMomentum( goodMomentumSeeds[i] );
00288 double pt = sqrt( (mom.x()*mom.x()) + (mom.y()*mom.y()) );
00289 if ( pt < 6. || pt > 2000. ) continue;
00290
00291
00292 seeds.push_back( goodMomentumSeeds[i] );
00293 findgoodMomentum = true;
00294 }
00295 if ( seeds.size() == 0 ) seeds = goodMomentumSeeds;
00296
00297 return findgoodMomentum;
00298 }
00299
00300
00301
00302 SeedContainer MuonSeedCleaner::SeedCandidates( std::vector<TrajectorySeed>& seeds, bool good ) {
00303
00304 SeedContainer theCandidate;
00305 theCandidate.clear();
00306
00307 bool longSeed = false;
00308 bool withFirstLayer = false ;
00309
00310
00311 for ( size_t i = 0; i < seeds.size(); i++ ) {
00312
00313 if (seeds[i].nHits() > 1 ) longSeed = true ;
00314
00315
00316 int idx = 0;
00317 for (edm::OwnVector<TrackingRecHit>::const_iterator r1 = seeds[i].recHits().first; r1 != seeds[i].recHits().second; r1++){
00318
00319 idx++;
00320 const GeomDet* gdet = theService->trackingGeometry()->idToDet( (*r1).geographicalId() );
00321 DetId geoId = gdet->geographicalId();
00322
00323 if ( geoId.subdetId() == MuonSubdetId::DT ) {
00324 DTChamberId DT_Id( (*r1).geographicalId() );
00325
00326 if (DT_Id.station() != 1) continue;
00327 withFirstLayer = true;
00328 }
00329 if ( geoId.subdetId() == MuonSubdetId::CSC ) {
00330 idx++;
00331 CSCDetId CSC_Id = CSCDetId( (*r1).geographicalId() );
00332
00333 if (CSC_Id.station() != 1) continue;
00334 withFirstLayer = true;
00335 }
00336 }
00337 bool goodseed = (longSeed && withFirstLayer) ? true : false ;
00338
00339 if ( goodseed == good ) theCandidate.push_back( seeds[i] );
00340 }
00341 return theCandidate;
00342
00343 }
00344
00345 std::vector<SeedContainer> MuonSeedCleaner::GroupSeeds( std::vector<TrajectorySeed>& seeds) {
00346
00347 std::vector<SeedContainer> seedCollection;
00348 seedCollection.clear();
00349 std::vector<TrajectorySeed> theGroup ;
00350 std::vector<bool> usedSeed(seeds.size(),false);
00351
00352
00353 for (size_t i= 0; i<seeds.size(); i++){
00354
00355 if (usedSeed[i]) continue;
00356 theGroup.push_back( seeds[i] );
00357 usedSeed[i] = true ;
00358
00359 GlobalPoint pos1 = SeedPosition( seeds[i] );
00360
00361 for (size_t j= i+1; j<seeds.size(); j++){
00362
00363
00364 unsigned int overlapping = OverlapSegments(seeds[i], seeds[j]) ;
00365 if ( !usedSeed[j] && overlapping > 0 ) {
00366
00367 if ( seeds[i].nHits() == overlapping && seeds[j].nHits() == overlapping ) {
00368 usedSeed[j] = true ;
00369 continue;
00370 }
00371 theGroup.push_back( seeds[j] );
00372 usedSeed[j] = true ;
00373 }
00374 if (usedSeed[j]) continue;
00375
00376
00377 GlobalPoint pos2 = SeedPosition( seeds[j]);
00378 double dh = pos1.eta() - pos2.eta() ;
00379 double df = pos1.phi() - pos2.phi() ;
00380 double dR = sqrt( (dh*dh) + (df*df) );
00381
00382 if ( dR > 0.3 && seeds[j].nHits() == 1) continue;
00383 if ( dR > 0.2 && seeds[j].nHits() > 1) continue;
00384 theGroup.push_back( seeds[j] );
00385 usedSeed[j] = true ;
00386 }
00387 sort(theGroup.begin(), theGroup.end(), lengthSorting ) ;
00388 seedCollection.push_back(theGroup);
00389
00390 theGroup.clear();
00391 }
00392 return seedCollection;
00393
00394 }
00395
00396 unsigned int MuonSeedCleaner::OverlapSegments( TrajectorySeed seed1, TrajectorySeed seed2 ) {
00397
00398 unsigned int overlapping = 0;
00399 for (edm::OwnVector<TrackingRecHit>::const_iterator r1 = seed1.recHits().first; r1 != seed1.recHits().second; r1++){
00400
00401 DetId id1 = (*r1).geographicalId();
00402 const GeomDet* gdet1 = theService->trackingGeometry()->idToDet( id1 );
00403 GlobalPoint gp1 = gdet1->toGlobal( (*r1).localPosition() );
00404
00405 for (edm::OwnVector<TrackingRecHit>::const_iterator r2 = seed2.recHits().first; r2 != seed2.recHits().second; r2++){
00406
00407 DetId id2 = (*r2).geographicalId();
00408 if (id1 != id2 ) continue;
00409
00410 const GeomDet* gdet2 = theService->trackingGeometry()->idToDet( id2 );
00411 GlobalPoint gp2 = gdet2->toGlobal( (*r2).localPosition() );
00412
00413 double dx = gp1.x() - gp2.x() ;
00414 double dy = gp1.y() - gp2.y() ;
00415 double dz = gp1.z() - gp2.z() ;
00416 double dL = sqrt( dx*dx + dy*dy + dz*dz);
00417
00418 if ( dL < 1. ) overlapping ++;
00419
00420 }
00421 }
00422 return overlapping ;
00423
00424 }
00425
00426 double MuonSeedCleaner::SeedChi2( TrajectorySeed seed ) {
00427
00428 double theChi2 = 0.;
00429 for (edm::OwnVector<TrackingRecHit>::const_iterator r1 = seed.recHits().first; r1 != seed.recHits().second; r1++){
00430
00431 theChi2 += NChi2OfSegment( *r1 );
00432 }
00433 theChi2 = theChi2 / seed.nHits() ;
00434
00435
00436 return theChi2 ;
00437
00438 }
00439
00440 int MuonSeedCleaner::SeedLength( TrajectorySeed seed ) {
00441
00442 int theHits = 0;
00443 for (edm::OwnVector<TrackingRecHit>::const_iterator r1 = seed.recHits().first; r1 != seed.recHits().second; r1++){
00444
00445 theHits += NRecHitsFromSegment( *r1 );
00446 }
00447
00448
00449 return theHits ;
00450
00451 }
00452
00453 GlobalPoint MuonSeedCleaner::SeedPosition( TrajectorySeed seed ) {
00454
00455
00456
00457 PTrajectoryStateOnDet pTSOD = seed.startingState();
00458 DetId SeedDetId(pTSOD.detId());
00459 const GeomDet* geoDet = theService->trackingGeometry()->idToDet( SeedDetId );
00460 TrajectoryStateOnSurface SeedTSOS = trajectoryStateTransform::transientState(pTSOD, &(geoDet->surface()), &*theService->magneticField());
00461 GlobalPoint pos = SeedTSOS.globalPosition();
00462
00463 return pos ;
00464
00465 }
00466
00467 GlobalVector MuonSeedCleaner::SeedMomentum( TrajectorySeed seed ) {
00468
00469
00470
00471 PTrajectoryStateOnDet pTSOD = seed.startingState();
00472 DetId SeedDetId(pTSOD.detId());
00473 const GeomDet* geoDet = theService->trackingGeometry()->idToDet( SeedDetId );
00474 TrajectoryStateOnSurface SeedTSOS = trajectoryStateTransform::transientState(pTSOD, &(geoDet->surface()), &*theService->magneticField());
00475 GlobalVector mom = SeedTSOS.globalMomentum();
00476
00477 return mom ;
00478
00479 }
00480
00481 int MuonSeedCleaner::NRecHitsFromSegment( const TrackingRecHit& rhit ) {
00482
00483 int NRechits = 0 ;
00484 const GeomDet* gdet = theService->trackingGeometry()->idToDet( rhit.geographicalId() );
00485 MuonTransientTrackingRecHit::MuonRecHitPointer theSeg =
00486 MuonTransientTrackingRecHit::specificBuild(gdet, rhit.clone() );
00487
00488 DetId geoId = gdet->geographicalId();
00489 if ( geoId.subdetId() == MuonSubdetId::DT ) {
00490 DTChamberId DT_Id( rhit.geographicalId() );
00491 std::vector<TrackingRecHit*> DThits = theSeg->recHits();
00492 int dt1DHits = 0;
00493 for (size_t j=0; j< DThits.size(); j++) {
00494 dt1DHits += (DThits[j]->recHits()).size();
00495 }
00496 NRechits = dt1DHits ;
00497 }
00498
00499 if ( geoId.subdetId() == MuonSubdetId::CSC ) {
00500 NRechits = (theSeg->recHits()).size() ;
00501 }
00502 return NRechits ;
00503 }
00504
00505 int MuonSeedCleaner::NRecHitsFromSegment( MuonTransientTrackingRecHit *rhit ) {
00506
00507 int NRechits = 0 ;
00508 DetId geoId = rhit->geographicalId();
00509 if ( geoId.subdetId() == MuonSubdetId::DT ) {
00510 DTChamberId DT_Id( geoId );
00511 std::vector<TrackingRecHit*> DThits = rhit->recHits();
00512 int dt1DHits = 0;
00513 for (size_t j=0; j< DThits.size(); j++) {
00514 dt1DHits += (DThits[j]->recHits()).size();
00515 }
00516 NRechits = dt1DHits ;
00517
00518 }
00519 if ( geoId.subdetId() == MuonSubdetId::CSC ) {
00520 NRechits = (rhit->recHits()).size() ;
00521
00522 }
00523 return NRechits;
00524
00525 }
00526
00527 double MuonSeedCleaner::NChi2OfSegment( const TrackingRecHit& rhit ) {
00528
00529 double NChi2 = 999999. ;
00530 const GeomDet* gdet = theService->trackingGeometry()->idToDet( rhit.geographicalId() );
00531 MuonTransientTrackingRecHit::MuonRecHitPointer theSeg =
00532 MuonTransientTrackingRecHit::specificBuild(gdet, rhit.clone() );
00533
00534 double dof = static_cast<double>( theSeg->degreesOfFreedom() );
00535 NChi2 = theSeg->chi2() / dof ;
00536
00537
00538 return NChi2 ;
00539 }
00540
00541