00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044 #include "DPGAnalysis/Skims/interface/CSCSkim.h"
00045
00046 using namespace std;
00047 using namespace edm;
00048
00049
00050
00051
00052
00053 CSCSkim::CSCSkim(const edm::ParameterSet& pset)
00054 {
00055
00056
00057 cscRecHitTag = pset.getParameter<edm::InputTag>("cscRecHitTag");
00058 cscSegmentTag = pset.getParameter<edm::InputTag>("cscSegmentTag");
00059 SAMuonTag = pset.getParameter<edm::InputTag>("SAMuonTag");
00060 GLBMuonTag = pset.getParameter<edm::InputTag>("GLBMuonTag");
00061 trackTag = pset.getParameter<edm::InputTag>("trackTag");
00062
00063
00064 bool isSimulation = false;
00065 isSimulation = pset.getUntrackedParameter<bool>("isSimulation",false);
00066 outputFileName = pset.getUntrackedParameter<std::string>("outputFileName","outputSkim.root");
00067 histogramFileName = pset.getUntrackedParameter<std::string>("histogramFileName","histos.root");
00068 typeOfSkim = pset.getUntrackedParameter<int>("typeOfSkim",1);
00069 nLayersWithHitsMinimum = pset.getUntrackedParameter<int>("nLayersWithHitsMinimum",3);
00070 minimumHitChambers = pset.getUntrackedParameter<int>("minimumHitChambers",1);
00071 minimumSegments = pset.getUntrackedParameter<int>("minimumSegments",3);
00072 demandChambersBothSides = pset.getUntrackedParameter<bool>("demandChambersBothSides",false);
00073 makeHistograms = pset.getUntrackedParameter<bool>("makeHistograms",false);
00074 makeHistogramsForMessyEvents = pset.getUntrackedParameter<bool>("makeHistogramsForMessyEvebts",false);
00075 whichEndcap = pset.getUntrackedParameter<int>("whichEndcap",2);
00076 whichStation = pset.getUntrackedParameter<int>("whichStation",3);
00077 whichRing = pset.getUntrackedParameter<int>("whichRing",2);
00078 whichChamber = pset.getUntrackedParameter<int>("whichChamber",24);
00079
00080
00081 pMin = pset.getUntrackedParameter<double>("pMin",3.);
00082 zLengthMin = pset.getUntrackedParameter<double>("zLengthMin",200.);
00083 nCSCHitsMin = pset.getUntrackedParameter<int>("nCSCHitsMin",9);
00084 zInnerMax = pset.getUntrackedParameter<double>("zInnerMax",9000.);
00085 nTrHitsMin = pset.getUntrackedParameter<int>("nTrHitsMin",8);
00086 zLengthTrMin = pset.getUntrackedParameter<double>("zLengthTrMin",180.);
00087 rExtMax = pset.getUntrackedParameter<double>("rExtMax",3000.);
00088 redChiSqMax = pset.getUntrackedParameter<double>("redChiSqMax",20.);
00089 nValidHitsMin = pset.getUntrackedParameter<int>("nValidHitsMin",8);
00090
00091 LogInfo("[CSCSkim] Setup")
00092 << "\n\t===== CSCSkim =====\n"
00093 << "\t\ttype of skim ...............................\t" << typeOfSkim
00094 << "\t\tminimum number of layers with hits .........\t" << nLayersWithHitsMinimum
00095 << "\n\t\tminimum number of chambers w/ hit layers..\t" << minimumHitChambers
00096 << "\n\t\tminimum number of segments ...............\t" << minimumSegments
00097 << "\n\t\tdemand chambers on both sides.............\t" << demandChambersBothSides
00098 << "\n\t\tmake histograms...........................\t" << makeHistograms
00099 << "\n\t\t..for messy events........................\t" << makeHistogramsForMessyEvents
00100 << "\n\t===================\n\n";
00101
00102 }
00103
00104
00105
00106
00107 CSCSkim::~CSCSkim()
00108 {
00109 }
00110
00111
00112
00113
00114
00115 void
00116 CSCSkim::beginJob()
00117 {
00118
00119 nEventsAnalyzed = 0;
00120 nEventsSelected = 0;
00121 nEventsChambersBothSides = 0;
00122 nEventsOverlappingChambers = 0;
00123 nEventsMessy = 0;
00124 nEventsCertainChamber = 0;
00125 nEventsDTOverlap = 0;
00126 nEventsHaloLike = 0;
00127 nEventsLongSATrack = 0;
00128 nEventsForBFieldStudies = 0;
00129 iRun = 0;
00130 iEvent = 0;
00131
00132 if (makeHistograms || makeHistogramsForMessyEvents) {
00133
00134
00135 theHistogramFile = new TFile(histogramFileName.c_str(), "RECREATE");
00136 theHistogramFile->cd();
00137
00138 if (makeHistograms) {
00139
00140 hxnRecHits = new TH1F("hxnRecHits","n RecHits",61,-0.5,60.5);
00141 hxnSegments = new TH1F("hxnSegments","n Segments",11,-0.5,10.5);
00142 hxnHitChambers = new TH1F("hxnHitsChambers","n chambers with hits",11,-0.5,10.5);
00143 hxnRecHitsSel = new TH1F("hxnRecHitsSel","n RecHits selected",61,-0.5,60.5);
00144
00145 xxP = new TH1F("xxP","P global",100,0.,200.);
00146 xxnValidHits = new TH1F("xxnValidHits","n valid hits global",61,-0.5,60.5);
00147 xxnTrackerHits = new TH1F("xxnTrackerHits","n tracker hits global",61,-0.5,60.5);
00148 xxnCSCHits = new TH1F("xxnCSCHits","n CSC hits global",41,-0.5,40.5);
00149 xxredChiSq = new TH1F("xxredChiSq","red chisq global",100,0.,100.);
00150
00151
00152
00153 }
00154 if (makeHistogramsForMessyEvents) {
00155
00156 mevnRecHits0 = new TH1F("mevnRecHits0","n RecHits",121,-0.5,120.5);
00157 mevnChambers0 = new TH1F("mevnChambers0","n chambers with hits",21,-0.5,20.5);
00158 mevnSegments0 = new TH1F("mevnSegments0","n Segments",21,-0.5,20.5);
00159 mevnRecHits1 = new TH1F("mevnRecHits1","n RecHits",100,0.,300.);
00160 mevnChambers1 = new TH1F("mevnChambers1","n chambers with hits",50,0.,50.);
00161 mevnSegments1 = new TH1F("mevnSegments1","n Segments",30,0.,30.);
00162 }
00163
00164 }
00165 }
00166
00167
00168
00169
00170 void
00171 CSCSkim::endJob() {
00172
00173
00174
00175 float fraction = 0.;
00176 if (nEventsAnalyzed > 0) {fraction = (float)nEventsSelected / (float)nEventsAnalyzed;}
00177
00178 LogInfo("[CSCSkim] Summary")
00179 << "\n\n\t====== CSCSkim ==========================================================\n"
00180 << "\t\ttype of skim ...............................\t" << typeOfSkim << "\n"
00181 << "\t\tevents analyzed ..............\t" << nEventsAnalyzed << "\n"
00182 << "\t\tevents selected ..............\t" << nEventsSelected << "\tfraction= " << fraction << std::endl
00183 << "\t\tevents chambers both sides ...\t" << nEventsChambersBothSides << "\n"
00184 << "\t\tevents w/ overlaps .......... \t" << nEventsOverlappingChambers << "\n"
00185 << "\t\tevents lots of hit chambers . \t" << nEventsMessy << "\n"
00186 << "\t\tevents from certain chamber . \t" << nEventsCertainChamber << "\n"
00187 << "\t\tevents in DT-CSC overlap .... \t" << nEventsDTOverlap << "\n"
00188 << "\t\tevents halo-like ............ \t" << nEventsHaloLike << "\n"
00189 << "\t\tevents w/ long SA track ..... \t" << nEventsLongSATrack << "\n"
00190 << "\t\tevents good for BField ..... \t" << nEventsForBFieldStudies << "\n"
00191 << "\t=========================================================================\n\n";
00192
00193 if (makeHistograms || makeHistogramsForMessyEvents) {
00194
00195 LogDebug("[CSCSkim]") << "======= write out my histograms ====\n" ;
00196 theHistogramFile->cd();
00197 if (makeHistograms) {
00198 hxnRecHits->Write();
00199 hxnSegments->Write();
00200 hxnHitChambers->Write();
00201 hxnRecHitsSel->Write();
00202 }
00203 if (makeHistogramsForMessyEvents) {
00204 mevnRecHits0->Write();
00205 mevnChambers0->Write();
00206 mevnSegments0->Write();
00207 mevnRecHits1->Write();
00208 mevnChambers1->Write();
00209 mevnSegments1->Write();
00210 }
00211 theHistogramFile->Close();
00212 }
00213 }
00214
00215
00216
00217
00218 bool
00219 CSCSkim::filter(edm::Event& event, const edm::EventSetup& eventSetup)
00220 {
00221
00222 nEventsAnalyzed++;
00223
00224 iRun = event.id().run();
00225 iEvent = event.id().event();
00226
00227 LogDebug("[CSCSkim] EventInfo") << "Run: " << iRun << "\tEvent: " << iEvent << "\tn Analyzed: " << nEventsAnalyzed;
00228
00229
00230 ESHandle<CSCGeometry> cscGeom;
00231 eventSetup.get<MuonGeometryRecord>().get(cscGeom);
00232
00233
00234 edm::Handle<CSCWireDigiCollection> wires;
00235 edm::Handle<CSCStripDigiCollection> strips;
00236
00237 if (isSimulation){
00238 event.getByLabel("simMuonCSCDigis","MuonCSCWireDigi",wires);
00239 event.getByLabel("simMuonCSCDigis","MuonCSCStripDigi",strips);
00240 }
00241 else {
00242 event.getByLabel("muonCSCDigis","MuonCSCWireDigi",wires);
00243 event.getByLabel("muonCSCDigis","MuonCSCStripDigi",strips);
00244 }
00245
00246
00247 Handle<CSCRecHit2DCollection> cscRecHits;
00248 event.getByLabel(cscRecHitTag,cscRecHits);
00249
00250
00251 Handle<CSCSegmentCollection> cscSegments;
00252 event.getByLabel(cscSegmentTag, cscSegments);
00253
00254
00255 Handle<reco::TrackCollection> saMuons;
00256 if (typeOfSkim == 8) {
00257 event.getByLabel(SAMuonTag,saMuons);
00258 }
00259
00260
00261 Handle<reco::TrackCollection> tracks;
00262 Handle<reco::MuonCollection> gMuons;
00263 if (typeOfSkim == 9) {
00264 event.getByLabel(SAMuonTag,saMuons);
00265 event.getByLabel(trackTag,tracks);
00266 event.getByLabel(GLBMuonTag,gMuons);
00267 }
00268
00269
00270
00271
00272
00273
00274
00275
00276 bool basicEvent = false;
00277 if (typeOfSkim == 1 || typeOfSkim == 2) {
00278 basicEvent = doCSCSkimming(cscRecHits,cscSegments);
00279 }
00280
00281
00282 bool goodOverlapEvent = false;
00283 if (typeOfSkim == 3) {
00284 goodOverlapEvent = doOverlapSkimming(cscSegments);
00285 if (goodOverlapEvent) {nEventsOverlappingChambers++;}
00286 }
00287
00288
00289 bool messyEvent = false;
00290 if (typeOfSkim == 4) {
00291 messyEvent = doMessyEventSkimming(cscRecHits,cscSegments);
00292 if (messyEvent) {nEventsMessy++;}
00293 }
00294
00295
00296 bool hasChamber = false;
00297 if (typeOfSkim == 5) {
00298 hasChamber = doCertainChamberSelection(wires,strips);
00299 if (hasChamber) {nEventsCertainChamber++;}
00300 }
00301
00302
00303 bool DTOverlapCandidate = false;
00304 if (typeOfSkim == 6) {
00305 DTOverlapCandidate = doDTOverlap(cscSegments);
00306 if (DTOverlapCandidate) {nEventsDTOverlap++;}
00307 }
00308
00309
00310 bool HaloLike = false;
00311 if (typeOfSkim == 7) {
00312 HaloLike = doHaloLike(cscSegments);
00313 if (HaloLike) {nEventsHaloLike++;}
00314 }
00315
00316
00317 bool LongSATrack = false;
00318 if (typeOfSkim == 8) {
00319 LongSATrack = doLongSATrack(saMuons);
00320 if (LongSATrack) {nEventsLongSATrack++;}
00321 }
00322
00323
00324 bool GoodForBFieldStudy = false;
00325 if (typeOfSkim == 9) {
00326 GoodForBFieldStudy = doBFieldStudySelection(saMuons,tracks,gMuons);
00327 if (GoodForBFieldStudy) {nEventsForBFieldStudies++;}
00328 }
00329
00330
00331
00332 bool selectThisEvent = false;
00333 if (typeOfSkim == 1 || typeOfSkim == 2) {selectThisEvent = basicEvent;}
00334 if (typeOfSkim == 3) {selectThisEvent = goodOverlapEvent;}
00335 if (typeOfSkim == 4) {selectThisEvent = messyEvent;}
00336 if (typeOfSkim == 5) {selectThisEvent = hasChamber;}
00337 if (typeOfSkim == 6) {selectThisEvent = DTOverlapCandidate;}
00338 if (typeOfSkim == 7) {selectThisEvent = HaloLike;}
00339 if (typeOfSkim == 8) {selectThisEvent = LongSATrack;}
00340 if (typeOfSkim == 9) {selectThisEvent = GoodForBFieldStudy;}
00341
00342 if (selectThisEvent) {nEventsSelected++;}
00343
00344 return selectThisEvent;
00345 }
00346
00347
00348
00349
00350
00351
00352
00353
00354 bool CSCSkim::doCSCSkimming(edm::Handle<CSCRecHit2DCollection> cscRecHits, edm::Handle<CSCSegmentCollection> cscSegments){
00355
00356
00357 int nRecHits = cscRecHits->size();
00358
00359
00360 int cntRecHit[600];
00361 for (int i = 0; i < 600; i++) {
00362 cntRecHit[i] = 0;
00363 }
00364
00365
00366
00367
00368
00369 CSCRecHit2DCollection::const_iterator recIt;
00370 for (recIt = cscRecHits->begin(); recIt != cscRecHits->end(); recIt++) {
00371
00372
00373 CSCDetId idrec = (CSCDetId)(*recIt).cscDetId();
00374 int kEndcap = idrec.endcap();
00375 int kRing = idrec.ring();
00376 int kStation = idrec.station();
00377 int kChamber = idrec.chamber();
00378 int kLayer = idrec.layer();
00379
00380
00381 int kSerial = chamberSerial( kEndcap, kStation, kRing, kChamber ) ;
00382
00383
00384
00385 int kDigit = (int) pow((float)10.,(float)(kLayer-1));
00386 cntRecHit[kSerial] += kDigit;
00387
00388 }
00389
00390
00391
00392
00393
00394
00395 int nChambersWithMinimalHits = 0;
00396 int nChambersWithMinimalHitsPOS = 0;
00397 int nChambersWithMinimalHitsNEG = 0;
00398 if (nRecHits > 0) {
00399 for (int i = 0; i < 600; i++) {
00400 if (cntRecHit[i] > 0) {
00401 int nLayersWithHits = 0;
00402 float dummy = (float) cntRecHit[i];
00403 for (int j = 5; j > -1; j--) {
00404 float digit = dummy / pow( (float)10., (float)j );
00405 int kCount = (int) digit;
00406 if (kCount > 0) nLayersWithHits++;
00407 dummy = dummy - ( (float) kCount) * pow( (float)10., (float)j );
00408 }
00409 if (nLayersWithHits > nLayersWithHitsMinimum) {
00410 if (i < 300) {nChambersWithMinimalHitsPOS++;}
00411 else {nChambersWithMinimalHitsNEG++;}
00412 }
00413 }
00414 }
00415 nChambersWithMinimalHits = nChambersWithMinimalHitsPOS + nChambersWithMinimalHitsNEG;
00416 }
00417
00418
00419 int nSegments = cscSegments->size();
00420
00421
00422
00423
00424
00425 if (makeHistograms) {
00426 hxnRecHits->Fill(nRecHits);
00427 if (nRecHits > 0) {
00428 hxnSegments->Fill(nSegments);
00429 hxnHitChambers->Fill(nChambersWithMinimalHits);
00430 }
00431 if (nChambersWithMinimalHits > 0) {
00432 hxnRecHitsSel->Fill(nRecHits);
00433 }
00434 }
00435
00436
00437
00438
00439 bool basicEvent = ( nChambersWithMinimalHits >= minimumHitChambers ) && ( nSegments >= minimumSegments );
00440
00441 bool chambersOnBothSides = ((nChambersWithMinimalHitsPOS >= minimumHitChambers) && (nChambersWithMinimalHitsNEG >= minimumHitChambers));
00442
00443 if (chambersOnBothSides) {nEventsChambersBothSides++;}
00444
00445 bool selectEvent = false;
00446 if (typeOfSkim == 1) {selectEvent = basicEvent;}
00447 if (typeOfSkim == 2) {selectEvent = chambersOnBothSides;}
00448
00449
00450
00451 LogDebug("[CSCSkim]") << "----- nRecHits = " << nRecHits
00452 << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits
00453 << "\tnSegments = " << nSegments
00454 << "\tselect? " << selectEvent << std::endl;
00455
00456
00457
00458
00459
00460
00461
00462
00463
00464
00465
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475 return selectEvent;
00476 }
00477
00478
00479
00480
00481
00482 bool CSCSkim::doOverlapSkimming(edm::Handle<CSCSegmentCollection> cscSegments){
00483
00484 const int nhitsMinimum = 4;
00485 const float chisqMaximum = 100.;
00486 const int nAllMaximum = 3;
00487
00488
00489
00490
00491
00492 int nAll[600];
00493 int nGood[600];
00494 for (int i=0; i<600; i++) {
00495 nAll[i] = 0;
00496 nGood[i] = 0;
00497 }
00498
00499
00500
00501
00502 for(CSCSegmentCollection::const_iterator it=cscSegments->begin(); it != cscSegments->end(); it++) {
00503
00504
00505 CSCDetId id = (CSCDetId)(*it).cscDetId();
00506 int kEndcap = id.endcap();
00507 int kStation = id.station();
00508 int kRing = id.ring();
00509 int kChamber = id.chamber();
00510 int kSerial = chamberSerial( kEndcap, kStation, kRing, kChamber);
00511
00512
00513 float chisq = (*it).chi2();
00514 int nhits = (*it).nRecHits();
00515
00516
00517 bool goodSegment = (nhits >= nhitsMinimum) && (chisq < chisqMaximum) ;
00518
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530 nAll[kSerial-1]++;
00531 if (goodSegment) nGood[kSerial]++;
00532
00533 }
00534
00535
00536
00537
00538
00539
00540 bool messyChamber = false;
00541 for (int i = 0; i < 600; i++) {
00542 if (nAll[i] > nAllMaximum) messyChamber = true;
00543 }
00544
00545
00546
00547 bool consecutiveChambers = false;
00548 for (int i = 0; i < 599; i++) {
00549 if ( (nGood[i]>0) && (nGood[i+1]>0) ) consecutiveChambers = true;
00550 }
00551
00552 bool selectThisEvent = !messyChamber && consecutiveChambers;
00553
00554 return selectThisEvent;
00555
00556 }
00557
00558
00559
00560
00561
00562
00563
00564 bool CSCSkim::doMessyEventSkimming(edm::Handle<CSCRecHit2DCollection> cscRecHits, edm::Handle<CSCSegmentCollection> cscSegments){
00565
00566
00567 int nRecHits = cscRecHits->size();
00568
00569
00570 int cntRecHit[600];
00571 for (int i = 0; i < 600; i++) {
00572 cntRecHit[i] = 0;
00573 }
00574
00575
00576
00577
00578
00579 CSCRecHit2DCollection::const_iterator recIt;
00580 for (recIt = cscRecHits->begin(); recIt != cscRecHits->end(); recIt++) {
00581
00582
00583 CSCDetId idrec = (CSCDetId)(*recIt).cscDetId();
00584 int kEndcap = idrec.endcap();
00585 int kRing = idrec.ring();
00586 int kStation = idrec.station();
00587 int kChamber = idrec.chamber();
00588 int kLayer = idrec.layer();
00589
00590
00591 int kSerial = chamberSerial( kEndcap, kStation, kRing, kChamber ) ;
00592
00593
00594
00595 int kDigit = (int) pow((float)10.,(float)(kLayer-1));
00596 cntRecHit[kSerial] += kDigit;
00597
00598 }
00599
00600
00601
00602
00603
00604
00605 int nChambersWithMinimalHits = 0;
00606 int nChambersWithMinimalHitsPOS = 0;
00607 int nChambersWithMinimalHitsNEG = 0;
00608 if (nRecHits > 0) {
00609 for (int i = 0; i < 600; i++) {
00610 if (cntRecHit[i] > 0) {
00611 int nLayersWithHits = 0;
00612 float dummy = (float) cntRecHit[i];
00613 for (int j = 5; j > -1; j--) {
00614 float digit = dummy / pow( (float)10., (float)j );
00615 int kCount = (int) digit;
00616 if (kCount > 0) nLayersWithHits++;
00617 dummy = dummy - ( (float) kCount) * pow( (float)10., (float)j );
00618 }
00619 if (nLayersWithHits > nLayersWithHitsMinimum) {
00620 if (i < 300) {nChambersWithMinimalHitsPOS++;}
00621 else {nChambersWithMinimalHitsNEG++;}
00622 }
00623 }
00624 }
00625 nChambersWithMinimalHits = nChambersWithMinimalHitsPOS + nChambersWithMinimalHitsNEG;
00626 }
00627
00628
00629 int nSegments = cscSegments->size();
00630
00631
00632
00633
00634
00635 if (makeHistogramsForMessyEvents) {
00636 if (nRecHits > 8) {
00637 mevnRecHits0->Fill(nRecHits);
00638 mevnChambers0->Fill(nChambersWithMinimalHits);
00639 mevnSegments0->Fill(nSegments);
00640 }
00641 if (nRecHits > 54) {
00642 double dummy = (double) nRecHits;
00643 if (dummy > 299.9) dummy = 299.9;
00644 mevnRecHits1->Fill(dummy);
00645 dummy = (double) nChambersWithMinimalHits;
00646 if (dummy > 49.9) dummy = 49.9;
00647 mevnChambers1->Fill(dummy);
00648 dummy = (double) nSegments;
00649 if (dummy > 29.9) dummy = 29.9;
00650 mevnSegments1->Fill(dummy);
00651 }
00652 }
00653
00654
00655
00656
00657
00658 bool selectEvent = false;
00659 if ( (nRecHits > 54) && (nChambersWithMinimalHits > 5) ) {selectEvent = true;}
00660
00661
00662 LogDebug("[CSCSkim]") << "----- nRecHits = " << nRecHits
00663 << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits
00664 << "\tnSegments = " << nSegments
00665 << "\tselect? " << selectEvent << std::endl;
00666
00667
00668
00669
00670
00671
00672
00673
00674
00675
00676
00677
00678
00679
00680
00681
00682
00683
00684
00685
00686 return selectEvent;
00687 }
00688
00689
00690
00691
00692
00693
00694
00695 bool CSCSkim::doCertainChamberSelection(edm::Handle<CSCWireDigiCollection> wires,
00696 edm::Handle<CSCStripDigiCollection> strips) {
00697
00698
00699 bool certainChamberIsPresentInWires = false;
00700 for (CSCWireDigiCollection::DigiRangeIterator jw=wires->begin(); jw!=wires->end(); jw++) {
00701 CSCDetId id = (CSCDetId)(*jw).first;
00702 int kEndcap = id.endcap();
00703 int kRing = id.ring();
00704 int kStation = id.station();
00705 int kChamber = id.chamber();
00706 if ( (kEndcap == whichEndcap) &&
00707 (kStation == whichStation) &&
00708 (kRing == whichRing) &&
00709 (kChamber == whichChamber) )
00710 {certainChamberIsPresentInWires = true;}
00711 }
00712
00713
00714
00715 bool certainChamberIsPresentInStrips = false;
00716 for (CSCStripDigiCollection::DigiRangeIterator js=strips->begin(); js!=strips->end(); js++) {
00717 CSCDetId id = (CSCDetId)(*js).first;
00718 int kEndcap = id.endcap();
00719 int kRing = id.ring();
00720 int kStation = id.station();
00721 int kChamber = id.chamber();
00722 if ( (kEndcap == whichEndcap) &&
00723 (kStation == whichStation) &&
00724 (kRing == whichRing) &&
00725 (kChamber == whichChamber) )
00726 {certainChamberIsPresentInStrips = true;}
00727 }
00728
00729 bool certainChamberIsPresent = certainChamberIsPresentInWires || certainChamberIsPresentInStrips;
00730
00731 return certainChamberIsPresent;
00732 }
00733
00734
00735
00736
00737
00738
00739
00740
00741 bool CSCSkim::doDTOverlap(Handle<CSCSegmentCollection> cscSegments) {
00742 const float chisqMax = 100.;
00743 const int nhitsMin = 5;
00744 const int maxNSegments = 3;
00745
00746
00747 bool DTOverlapCandidate = false;
00748 int cntMEP13[36];
00749 int cntMEN13[36];
00750 int cntMEP22[36];
00751 int cntMEN22[36];
00752 int cntMEP32[36];
00753 int cntMEN32[36];
00754 for (int i=0; i<36; ++i) {
00755 cntMEP13[i] = 0;
00756 cntMEN13[i] = 0;
00757 cntMEP22[i] = 0;
00758 cntMEN22[i] = 0;
00759 cntMEP32[i] = 0;
00760 cntMEN32[i] = 0;
00761 }
00762
00763
00764
00765
00766
00767 int nSegments = cscSegments->size();
00768 if (nSegments < 2) return DTOverlapCandidate;
00769
00770 for(CSCSegmentCollection::const_iterator it=cscSegments->begin(); it != cscSegments->end(); it++) {
00771
00772 CSCDetId id = (CSCDetId)(*it).cscDetId();
00773 int kEndcap = id.endcap();
00774 int kStation = id.station();
00775 int kRing = id.ring();
00776 int kChamber = id.chamber();
00777
00778 float chisq = (*it).chi2();
00779 int nhits = (*it).nRecHits();
00780 bool goodSegment = (chisq < chisqMax) && (nhits >= nhitsMin) ;
00781 if (goodSegment) {
00782 if ( (kStation == 1) && (kRing == 3) ) {
00783 if (kEndcap == 1) {cntMEP13[kChamber-1]++;}
00784 if (kEndcap == 2) {cntMEN13[kChamber-1]++;}
00785 }
00786 if ( (kStation == 2) && (kRing == 2) ) {
00787 if (kEndcap == 1) {cntMEP22[kChamber-1]++;}
00788 if (kEndcap == 2) {cntMEN22[kChamber-1]++;}
00789 }
00790 if ( (kStation == 3) && (kRing == 2) ) {
00791 if (kEndcap == 1) {cntMEP32[kChamber-1]++;}
00792 if (kEndcap == 2) {cntMEN32[kChamber-1]++;}
00793 }
00794 }
00795 }
00796
00797
00798
00799
00800 bool tooManySegments = false;
00801 for (int i=0; i<36; ++i) {
00802 if ( (cntMEP13[i] > maxNSegments) ||
00803 (cntMEN13[i] > maxNSegments) ||
00804 (cntMEP22[i] > maxNSegments) ||
00805 (cntMEN22[i] > maxNSegments) ||
00806 (cntMEP32[i] > maxNSegments) ||
00807 (cntMEN32[i] > maxNSegments) ) tooManySegments = true;
00808 }
00809 if (tooManySegments) {
00810 return DTOverlapCandidate;
00811 }
00812
00813
00814
00815
00816 bool matchup = false;
00817 for (int i=0; i<36; ++i) {
00818 if ( (cntMEP13[i] > 0) && (cntMEP22[i]+cntMEP32[i] > 0) ) {matchup = true;}
00819 if ( (cntMEN13[i] > 0) && (cntMEN22[i]+cntMEN32[i] > 0) ) {matchup = true;}
00820 }
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832
00833
00834
00835
00836 DTOverlapCandidate = matchup;
00837 return DTOverlapCandidate;
00838 }
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848
00849 bool CSCSkim::doHaloLike(Handle<CSCSegmentCollection> cscSegments) {
00850 const float chisqMax = 100.;
00851 const int nhitsMin = 5;
00852 const int maxNSegments = 3;
00853
00854
00855 bool HaloLike = false;
00856 int cntMEP11[36];
00857 int cntMEN11[36];
00858 int cntMEP12[36];
00859 int cntMEN12[36];
00860 int cntMEP21[36];
00861 int cntMEN21[36];
00862 int cntMEP31[36];
00863 int cntMEN31[36];
00864 int cntMEP41[36];
00865 int cntMEN41[36];
00866 for (int i=0; i<36; ++i) {
00867 cntMEP11[i] = 0;
00868 cntMEN11[i] = 0;
00869 cntMEP12[i] = 0;
00870 cntMEN12[i] = 0;
00871 cntMEP21[i] = 0;
00872 cntMEN21[i] = 0;
00873 cntMEP31[i] = 0;
00874 cntMEN31[i] = 0;
00875 cntMEP41[i] = 0;
00876 cntMEN41[i] = 0;
00877 }
00878
00879
00880
00881
00882 int nSegments = cscSegments->size();
00883 if (nSegments < 4) return HaloLike;
00884
00885 for(CSCSegmentCollection::const_iterator it=cscSegments->begin(); it != cscSegments->end(); it++) {
00886
00887 CSCDetId id = (CSCDetId)(*it).cscDetId();
00888 int kEndcap = id.endcap();
00889 int kStation = id.station();
00890 int kRing = id.ring();
00891 int kChamber = id.chamber();
00892
00893 float chisq = (*it).chi2();
00894 int nhits = (*it).nRecHits();
00895 bool goodSegment = (chisq < chisqMax) && (nhits >= nhitsMin) ;
00896 if (goodSegment) {
00897 if ( (kStation == 1) && (kRing == 1) ) {
00898 if (kEndcap == 1) {cntMEP11[kChamber-1]++;}
00899 if (kEndcap == 2) {cntMEN11[kChamber-1]++;}
00900 }
00901 if ( (kStation == 1) && (kRing == 2) ) {
00902 if (kEndcap == 1) {cntMEP12[kChamber-1]++;}
00903 if (kEndcap == 2) {cntMEN12[kChamber-1]++;}
00904 }
00905 if ( (kStation == 2) && (kRing == 1) ) {
00906 if (kEndcap == 1) {cntMEP21[kChamber-1]++;}
00907 if (kEndcap == 2) {cntMEN21[kChamber-1]++;}
00908 }
00909 if ( (kStation == 3) && (kRing == 1) ) {
00910 if (kEndcap == 1) {cntMEP31[kChamber-1]++;}
00911 if (kEndcap == 2) {cntMEN31[kChamber-1]++;}
00912 }
00913 if ( (kStation == 4) && (kRing == 1) ) {
00914 if (kEndcap == 1) {cntMEP41[kChamber-1]++;}
00915 if (kEndcap == 2) {cntMEN41[kChamber-1]++;}
00916 }
00917 }
00918 }
00919
00920
00921
00922
00923 bool tooManySegments = false;
00924 for (int i=0; i<36; ++i) {
00925 if ( (cntMEP11[i] > 3*maxNSegments) ||
00926 (cntMEN11[i] > 3*maxNSegments) ||
00927 (cntMEP12[i] > maxNSegments) ||
00928 (cntMEN12[i] > maxNSegments) ||
00929 (cntMEP21[i] > maxNSegments) ||
00930 (cntMEN21[i] > maxNSegments) ||
00931 (cntMEP31[i] > maxNSegments) ||
00932 (cntMEN31[i] > maxNSegments) ||
00933 (cntMEP41[i] > maxNSegments) ||
00934 (cntMEN41[i] > maxNSegments) ) tooManySegments = true;
00935 }
00936 if (tooManySegments) {
00937 return HaloLike;
00938 }
00939
00940
00941
00942
00943 bool matchup = false;
00944 for (int i=0; i<36; ++i) {
00945 if ( (cntMEP11[i]+cntMEP12[i] > 0) &&
00946 (cntMEP21[i] > 0) &&
00947 (cntMEP31[i] > 0) &&
00948 (cntMEP41[i] > 0) ) {matchup = true;}
00949 if ( (cntMEN11[i]+cntMEN12[i] > 0) &&
00950 (cntMEN21[i] > 0) &&
00951 (cntMEN31[i] > 0) &&
00952 (cntMEN41[i] > 0) ) {matchup = true;}
00953 }
00954
00955
00956
00957
00958
00959
00960
00961
00962
00963
00964
00965
00966
00967
00968
00969
00970
00971
00972
00973
00974
00975
00976
00977
00978
00979
00980
00981
00982
00983
00984
00985
00986 HaloLike = matchup;
00987 return HaloLike;
00988 }
00989
00990
00991
00992
00993
00994 bool CSCSkim::doLongSATrack(edm::Handle<reco::TrackCollection> saMuons) {
00995
00996 const float zDistanceMax = 2500.;
00997 const float zDistanceMin = 700.;
00998 const int nCSCHitsMin = 25;
00999 const int nCSCHitsMax = 50;
01000 const float zInnerMax = 80000.;
01001
01002 const int nNiceMuonsMin = 1;
01003
01004
01005
01006
01007
01008 int nNiceMuons = 0;
01009
01010 for(reco::TrackCollection::const_iterator muon = saMuons->begin(); muon != saMuons->end(); ++ muon ) {
01011
01012
01013 math::XYZVector innerMo = muon->innerMomentum();
01014 GlobalVector im(innerMo.x(),innerMo.y(),innerMo.z());
01015 math::XYZPoint innerPo = muon->innerPosition();
01016 GlobalPoint ip(innerPo.x(), innerPo.y(),innerPo.z());
01017 math::XYZPoint outerPo = muon->outerPosition();
01018 GlobalPoint op(outerPo.x(), outerPo.y(),outerPo.z());
01019 float zInner = ip.z();
01020 float zOuter = op.z();
01021 float zDistance = fabs(zOuter-zInner);
01022
01023
01024
01025
01026 int nDTHits = 0;
01027 int nCSCHits = 0;
01028 for (trackingRecHit_iterator hit = muon->recHitsBegin(); hit != muon->recHitsEnd(); ++hit ) {
01029 const DetId detId( (*hit)->geographicalId() );
01030 if (detId.det() == DetId::Muon) {
01031 if (detId.subdetId() == MuonSubdetId::DT) {
01032
01033
01034 nDTHits++;
01035 }
01036 else if (detId.subdetId() == MuonSubdetId::CSC) {
01037
01038
01039 nCSCHits++;
01040 }
01041 }
01042 }
01043
01044
01045 if ( (zDistance < zDistanceMax) && (zDistance > zDistanceMin)
01046 && (nCSCHits > nCSCHitsMin) && (nCSCHits < nCSCHitsMax)
01047 && (min ( fabs(zInner), fabs(zOuter) ) < zInnerMax)
01048 && (fabs(innerMo.z()) > 0.000000001) ) {
01049 nNiceMuons++;
01050 }
01051 }
01052
01053 bool select = (nNiceMuons >= nNiceMuonsMin);
01054
01055 return select;
01056 }
01057
01058
01059
01060
01061
01062
01063
01064
01065
01066
01067
01068
01069
01070 bool CSCSkim::doBFieldStudySelection(edm::Handle<reco::TrackCollection> saMuons, edm::Handle<reco::TrackCollection> tracks, edm::Handle<reco::MuonCollection> gMuons) {
01071
01072 bool acceptThisEvent = false;
01073
01074
01075
01076
01077 int nGoodSAMuons = 0;
01078 for (reco::TrackCollection::const_iterator muon = saMuons->begin(); muon != saMuons->end(); ++ muon ) {
01079 float preco = muon->p();
01080
01081 math::XYZPoint innerPo = muon->innerPosition();
01082 GlobalPoint iPnt(innerPo.x(), innerPo.y(),innerPo.z());
01083 math::XYZPoint outerPo = muon->outerPosition();
01084 GlobalPoint oPnt(outerPo.x(), outerPo.y(),outerPo.z());
01085 float zLength = abs( iPnt.z() - oPnt.z() );
01086
01087 math::XYZVector innerMom = muon->innerMomentum();
01088 GlobalVector iP(innerMom.x(), innerMom.y(), innerMom.z() );
01089 math::XYZVector outerMom = muon->outerMomentum();
01090 GlobalVector oP(outerMom.x(), outerMom.y(), outerMom.z() );
01091
01092 const float zRef = 300.;
01093 float xExt = 10000.;
01094 float yExt = 10000.;
01095 if (abs(oPnt.z()) < abs(iPnt.z())) {
01096 float deltaZ = 0.;
01097 if (oPnt.z() > 0) {
01098 deltaZ = zRef - oPnt.z();
01099 } else {
01100 deltaZ = -zRef - oPnt.z();
01101 }
01102 xExt = oPnt.x() + deltaZ * oP.x() / oP.z();
01103 yExt = oPnt.y() + deltaZ * oP.y() / oP.z();
01104 } else {
01105 float deltaZ = 0.;
01106 if (iPnt.z() > 0) {
01107 deltaZ = zRef - iPnt.z();
01108 } else {
01109 deltaZ = -zRef - iPnt.z();
01110 }
01111 xExt = iPnt.x() + deltaZ * iP.x() / iP.z();
01112 yExt = iPnt.y() + deltaZ * iP.y() / iP.z();
01113 }
01114 float rExt = sqrt( xExt*xExt + yExt*yExt );
01115
01116 int kHit = 0;
01117 int nDTHits = 0;
01118 int nCSCHits = 0;
01119 for (trackingRecHit_iterator hit = muon->recHitsBegin(); hit != muon->recHitsEnd(); ++hit ) {
01120 ++kHit;
01121 const DetId detId( (*hit)->geographicalId() );
01122 if (detId.det() == DetId::Muon) {
01123 if (detId.subdetId() == MuonSubdetId::DT) {
01124 nDTHits++;
01125 }
01126 else if (detId.subdetId() == MuonSubdetId::CSC) {
01127 nCSCHits++;
01128 }
01129 }
01130 }
01131
01132 float zInner = -1.;
01133 float zOuter = 20000.;
01134 if (nCSCHits >= nCSCHitsMin) {
01135 if (abs(iPnt.z()) < abs(iPnt.z())) {
01136 zInner = iPnt.z();
01137 zOuter = oPnt.z();
01138 } else {
01139 zInner = oPnt.z();
01140 zOuter = iPnt.z();
01141 }
01142 }
01143
01144 bool goodSAMuon = (preco > pMin)
01145 && ( zLength > zLengthMin )
01146 && ( nCSCHits >= nCSCHitsMin )
01147 && ( zInner < zInnerMax )
01148 && ( rExt < rExtMax ) ;
01149
01150 if (goodSAMuon) {nGoodSAMuons++;}
01151
01152 }
01153
01154
01155
01156
01157
01158
01159
01160 int nGoodTracks = 0;
01161 for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); ++ track ) {
01162 float preco = track->p();
01163 int n = track->recHitsSize();
01164
01165 math::XYZPoint innerPo = track->innerPosition();
01166 GlobalPoint iPnt(innerPo.x(), innerPo.y(),innerPo.z());
01167 math::XYZPoint outerPo = track->outerPosition();
01168 GlobalPoint oPnt(outerPo.x(), outerPo.y(),outerPo.z());
01169 float zLength = abs( iPnt.z() - oPnt.z() );
01170
01171 math::XYZVector innerMom = track->innerMomentum();
01172 GlobalVector iP(innerMom.x(), innerMom.y(), innerMom.z() );
01173 math::XYZVector outerMom = track->outerMomentum();
01174 GlobalVector oP(outerMom.x(), outerMom.y(), outerMom.z() );
01175
01176 const float zRef = 300.;
01177 float xExt = 10000.;
01178 float yExt = 10000.;
01179 if (abs(oPnt.z()) > abs(iPnt.z())) {
01180 float deltaZ = 0.;
01181 if (oPnt.z() > 0) {
01182 deltaZ = zRef - oPnt.z();
01183 } else {
01184 deltaZ = -zRef - oPnt.z();
01185 }
01186 xExt = oPnt.x() + deltaZ * oP.x() / oP.z();
01187 yExt = oPnt.y() + deltaZ * oP.y() / oP.z();
01188 } else {
01189 float deltaZ = 0.;
01190 if (iPnt.z() > 0) {
01191 deltaZ = zRef - iPnt.z();
01192 } else {
01193 deltaZ = -zRef - iPnt.z();
01194 }
01195 xExt = iPnt.x() + deltaZ * iP.x() / iP.z();
01196 yExt = iPnt.y() + deltaZ * iP.y() / iP.z();
01197 }
01198 float rExt = sqrt( xExt*xExt + yExt*yExt );
01199
01200 float zInner = -1.;
01201 float zOuter = 20000.;
01202 if (n >= nTrHitsMin) {
01203 if (abs(iPnt.z()) < abs(iPnt.z())) {
01204 zInner = iPnt.z();
01205 zOuter = oPnt.z();
01206 } else {
01207 zInner = oPnt.z();
01208 zOuter = iPnt.z();
01209 }
01210 }
01211
01212 bool goodTrack = (preco > pMin)
01213 && (n >= nTrHitsMin)
01214 && (zLength > zLengthTrMin)
01215 && ( rExt < rExtMax ) ;
01216
01217 if (goodTrack) {nGoodTracks++;}
01218
01219 }
01220
01221
01222
01223
01224
01225 int nGoodGlobalMuons = 0;
01226 for (reco::MuonCollection::const_iterator global = gMuons->begin(); global != gMuons->end(); ++global ) {
01227
01228 if (global->isGlobalMuon()) {
01229
01230 float pDef = global->p();
01231 float redChiSq = global->globalTrack()->normalizedChi2();
01232 const reco::HitPattern& hp = (global->globalTrack())->hitPattern();
01233
01234
01235 int nTrackerHits = hp.numberOfValidTrackerHits();
01236
01237
01238
01239 int nDTHits = 0;
01240 int nCSCHits = 0;
01241 for (trackingRecHit_iterator hit = (global->globalTrack())->recHitsBegin(); hit != (global->globalTrack())->recHitsEnd(); ++hit ) {
01242 const DetId detId( (*hit)->geographicalId() );
01243 if (detId.det() == DetId::Muon) {
01244 if (detId.subdetId() == MuonSubdetId::DT) {
01245 nDTHits++;
01246 }
01247 else if (detId.subdetId() == MuonSubdetId::CSC) {
01248 nCSCHits++;
01249 }
01250 }
01251 }
01252
01253 bool goodGlobalMuon = (pDef > pMin)
01254 && ( nTrackerHits >= nValidHitsMin )
01255 && ( nCSCHits >= nCSCHitsMin )
01256 && ( redChiSq < redChiSqMax );
01257
01258 if (goodGlobalMuon) {nGoodGlobalMuons++;}
01259
01260 }
01261 }
01262
01263
01264
01265
01266
01267 acceptThisEvent = ( (nGoodSAMuons > 0) && (nGoodTracks > 0) ) || (nGoodGlobalMuons > 0) ;
01268
01269 return acceptThisEvent;
01270 }
01271
01272
01273
01274
01275
01276 int CSCSkim::chamberSerial( int kEndcap, int kStation, int kRing, int kChamber ) {
01277 int kSerial = kChamber;
01278 if (kStation == 1 && kRing == 1) {kSerial = kChamber;}
01279 if (kStation == 1 && kRing == 2) {kSerial = kChamber + 36;}
01280 if (kStation == 1 && kRing == 3) {kSerial = kChamber + 72;}
01281 if (kStation == 1 && kRing == 4) {kSerial = kChamber;}
01282 if (kStation == 2 && kRing == 1) {kSerial = kChamber + 108;}
01283 if (kStation == 2 && kRing == 2) {kSerial = kChamber + 126;}
01284 if (kStation == 3 && kRing == 1) {kSerial = kChamber + 162;}
01285 if (kStation == 3 && kRing == 2) {kSerial = kChamber + 180;}
01286 if (kStation == 4 && kRing == 1) {kSerial = kChamber + 216;}
01287 if (kStation == 4 && kRing == 2) {kSerial = kChamber + 234;}
01288 if (kEndcap == 2) {kSerial = kSerial + 300;}
01289 return kSerial;
01290 }
01291
01292
01293 DEFINE_FWK_MODULE(CSCSkim);