CMS 3D CMS Logo

CSCSkim.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CSCSkim
4 // Class: CSCSkim
5 //
13 //
14 // Original Author: Michael Schmitt
15 // Created: Sat Jul 12 17:43:33 CEST 2008
16 //
17 //
18 //======================================================================
19 //
20 // CSCSkim:
21 //
22 // A simple skim module for extracting generally useful events from
23 // the cosmic-ray runs (CRUZET-n and CRAFT). The selected events
24 // should be useful for most CSC-DPG and muon-POG studies. However,
25 // the selection requirements may bias the sample with respect to
26 // trigger requirements and certain noise and efficiency-related
27 // studies.
28 //
29 // Types of skims: (typeOfSkim control word)
30 // 1 = loose skim demanding hit chambers and/or segments
31 // 2 = ask for hit chambers in both endcaps
32 // 3 = segments in neighboring chambers - good for alignment
33 // 4 = messy events
34 // 5 = select events with DIGIs from one particular chamber
35 // 6 = overlap with DT
36 // 7 = nearly horizontal track going through ME1/1,2/1,3/1,4/1
37 // 8 = ask for one long cosmic stand-alone muon track
38 // 9 = selection for magnetic field studies
39 //
40 //
41 //======================================================================
42 
44 
45 using namespace std;
46 using namespace edm;
47 
48 //===================
49 // CONSTRUCTOR
50 //===================
51 CSCSkim::CSCSkim(const edm::ParameterSet& pset) : m_CSCGeomToken(esConsumes()) {
52  // tokens from tags
53 
54  // Really should define the wire and digi tags in config, but for now, to avoid having to modify
55  // multiple config files, just hard-code those tags, to be equivalent to the pre-consumes code
56 
57  // wds_token = consumes<CSCWireDigiCollection>(pset.getParameter<InputTag>("simWireDigiTag"));
58  // sds_token = consumes<CSCStripDigiCollection>(pset.getParameter<InputTag>("simStripDigiTag"));
59  // wdr_token = consumes<CSCWireDigiCollection>(pset.getParameter<InputTag>("wireDigiTag"));
60  // sdr_token = consumes<CSCStripDigiCollection>(pset.getParameter<InputTag>("stripDigiTag"));
61 
62  wds_token = consumes<CSCWireDigiCollection>(edm::InputTag("simMuonCSCDigis", "MuonCSCWireDigi"));
63  sds_token = consumes<CSCStripDigiCollection>(edm::InputTag("simMuonCSCDigis", "MuonCSCStripDigi"));
64  wdr_token = consumes<CSCWireDigiCollection>(edm::InputTag("muonCSCDigis", "MuonCSCWireDigi"));
65  sdr_token = consumes<CSCStripDigiCollection>(edm::InputTag("muonCSCDigis", "MuonCSCStripDigi"));
66 
67  rh_token = consumes<CSCRecHit2DCollection>(pset.getParameter<InputTag>("cscRecHitTag"));
68  seg_token = consumes<CSCSegmentCollection>(pset.getParameter<InputTag>("cscSegmentTag"));
69 
70  sam_token = consumes<reco::TrackCollection>(pset.getParameter<InputTag>("SAMuonTag"));
71  trk_token = consumes<reco::TrackCollection>(pset.getParameter<InputTag>("trackTag"));
72  glm_token = consumes<reco::MuonCollection>(pset.getParameter<InputTag>("GLBMuonTag"));
73 
74  // Get the various input parameters
75  outputFileName = pset.getUntrackedParameter<std::string>("outputFileName", "outputSkim.root");
76  histogramFileName = pset.getUntrackedParameter<std::string>("histogramFileName", "histos.root");
77  typeOfSkim = pset.getUntrackedParameter<int>("typeOfSkim", 1);
78  nLayersWithHitsMinimum = pset.getUntrackedParameter<int>("nLayersWithHitsMinimum", 3);
79  minimumHitChambers = pset.getUntrackedParameter<int>("minimumHitChambers", 1);
80  minimumSegments = pset.getUntrackedParameter<int>("minimumSegments", 3);
81  demandChambersBothSides = pset.getUntrackedParameter<bool>("demandChambersBothSides", false);
82  makeHistograms = pset.getUntrackedParameter<bool>("makeHistograms", false);
83  makeHistogramsForMessyEvents = pset.getUntrackedParameter<bool>("makeHistogramsForMessyEvebts", false);
84  whichEndcap = pset.getUntrackedParameter<int>("whichEndcap", 2);
85  whichStation = pset.getUntrackedParameter<int>("whichStation", 3);
86  whichRing = pset.getUntrackedParameter<int>("whichRing", 2);
87  whichChamber = pset.getUntrackedParameter<int>("whichChamber", 24);
88 
89  // for BStudy selection (skim type 9)
90  pMin = pset.getUntrackedParameter<double>("pMin", 3.);
91  zLengthMin = pset.getUntrackedParameter<double>("zLengthMin", 200.);
92  nCSCHitsMin = pset.getUntrackedParameter<int>("nCSCHitsMin", 9);
93  zInnerMax = pset.getUntrackedParameter<double>("zInnerMax", 9000.);
94  nTrHitsMin = pset.getUntrackedParameter<int>("nTrHitsMin", 8);
95  zLengthTrMin = pset.getUntrackedParameter<double>("zLengthTrMin", 180.);
96  rExtMax = pset.getUntrackedParameter<double>("rExtMax", 3000.);
97  redChiSqMax = pset.getUntrackedParameter<double>("redChiSqMax", 20.);
98  nValidHitsMin = pset.getUntrackedParameter<int>("nValidHitsMin", 8);
99 
100  LogInfo("[CSCSkim] Setup") << "\n\t===== CSCSkim =====\n"
101  << "\t\ttype of skim ...............................\t" << typeOfSkim
102  << "\t\tminimum number of layers with hits .........\t" << nLayersWithHitsMinimum
103  << "\n\t\tminimum number of chambers w/ hit layers..\t" << minimumHitChambers
104  << "\n\t\tminimum number of segments ...............\t" << minimumSegments
105  << "\n\t\tdemand chambers on both sides.............\t" << demandChambersBothSides
106  << "\n\t\tmake histograms...........................\t" << makeHistograms
107  << "\n\t\t..for messy events........................\t" << makeHistogramsForMessyEvents
108  << "\n\t===================\n\n";
109 }
110 
111 //===================
112 // DESTRUCTOR
113 //===================
115 
116 //================
117 // BEGIN JOB
118 //================
120  // set counters to zero
121  nEventsAnalyzed = 0;
122  nEventsSelected = 0;
125  nEventsMessy = 0;
127  nEventsDTOverlap = 0;
128  nEventsHaloLike = 0;
129  nEventsLongSATrack = 0;
131  iRun = 0;
132  iEvent = 0;
133 
135  // Create the root file for the histograms
136  theHistogramFile = new TFile(histogramFileName.c_str(), "RECREATE");
137  theHistogramFile->cd();
138 
139  if (makeHistograms) {
140  // book histograms for the skimming module
141  hxnRecHits = new TH1F("hxnRecHits", "n RecHits", 61, -0.5, 60.5);
142  hxnSegments = new TH1F("hxnSegments", "n Segments", 11, -0.5, 10.5);
143  hxnHitChambers = new TH1F("hxnHitsChambers", "n chambers with hits", 11, -0.5, 10.5);
144  hxnRecHitsSel = new TH1F("hxnRecHitsSel", "n RecHits selected", 61, -0.5, 60.5);
145 
146  xxP = new TH1F("xxP", "P global", 100, 0., 200.);
147  xxnValidHits = new TH1F("xxnValidHits", "n valid hits global", 61, -0.5, 60.5);
148  xxnTrackerHits = new TH1F("xxnTrackerHits", "n tracker hits global", 61, -0.5, 60.5);
149  xxnCSCHits = new TH1F("xxnCSCHits", "n CSC hits global", 41, -0.5, 40.5);
150  xxredChiSq = new TH1F("xxredChiSq", "red chisq global", 100, 0., 100.);
151  }
153  // book histograms for the messy event skimming module
154  mevnRecHits0 = new TH1F("mevnRecHits0", "n RecHits", 121, -0.5, 120.5);
155  mevnChambers0 = new TH1F("mevnChambers0", "n chambers with hits", 21, -0.5, 20.5);
156  mevnSegments0 = new TH1F("mevnSegments0", "n Segments", 21, -0.5, 20.5);
157  mevnRecHits1 = new TH1F("mevnRecHits1", "n RecHits", 100, 0., 300.);
158  mevnChambers1 = new TH1F("mevnChambers1", "n chambers with hits", 50, 0., 50.);
159  mevnSegments1 = new TH1F("mevnSegments1", "n Segments", 30, 0., 30.);
160  }
161  }
162 }
163 
164 //================
165 // END JOB
166 //================
168  // Write out results
169 
170  float fraction = 0.;
171  if (nEventsAnalyzed > 0) {
173  }
174 
175  LogInfo("[CSCSkim] Summary") << "\n\n\t====== CSCSkim ==========================================================\n"
176  << "\t\ttype of skim ...............................\t" << typeOfSkim << "\n"
177  << "\t\tevents analyzed ..............\t" << nEventsAnalyzed << "\n"
178  << "\t\tevents selected ..............\t" << nEventsSelected
179  << "\tfraction= " << fraction << std::endl
180  << "\t\tevents chambers both sides ...\t" << nEventsChambersBothSides << "\n"
181  << "\t\tevents w/ overlaps .......... \t" << nEventsOverlappingChambers << "\n"
182  << "\t\tevents lots of hit chambers . \t" << nEventsMessy << "\n"
183  << "\t\tevents from certain chamber . \t" << nEventsCertainChamber << "\n"
184  << "\t\tevents in DT-CSC overlap .... \t" << nEventsDTOverlap << "\n"
185  << "\t\tevents halo-like ............ \t" << nEventsHaloLike << "\n"
186  << "\t\tevents w/ long SA track ..... \t" << nEventsLongSATrack << "\n"
187  << "\t\tevents good for BField ..... \t" << nEventsForBFieldStudies << "\n"
188  << "\t=========================================================================\n\n";
189 
191  // Write the histos to file
192  LogDebug("[CSCSkim]") << "======= write out my histograms ====\n";
193  theHistogramFile->cd();
194  if (makeHistograms) {
195  hxnRecHits->Write();
196  hxnSegments->Write();
197  hxnHitChambers->Write();
198  hxnRecHitsSel->Write();
199  }
201  mevnRecHits0->Write();
202  mevnChambers0->Write();
203  mevnSegments0->Write();
204  mevnRecHits1->Write();
205  mevnChambers1->Write();
206  mevnSegments1->Write();
207  }
208  theHistogramFile->Close();
209  }
210 }
211 
212 //================
213 // FILTER MAIN
214 //================
216  // increment counter
217  nEventsAnalyzed++;
218 
219  iRun = event.id().run();
220  iEvent = event.id().event();
221 
222  LogDebug("[CSCSkim] EventInfo") << "Run: " << iRun << "\tEvent: " << iEvent << "\tn Analyzed: " << nEventsAnalyzed;
223 
224  // Get the CSC Geometry :
225  ESHandle<CSCGeometry> cscGeom = eventSetup.getHandle(m_CSCGeomToken);
226 
227  // Get the DIGI collections
230 
231  if (event.eventAuxiliary().isRealData()) {
232  event.getByToken(wdr_token, wires);
233  event.getByToken(sdr_token, strips);
234  } else {
235  event.getByToken(wds_token, wires);
236  event.getByToken(sds_token, strips);
237  }
238 
239  // Get the RecHits collection :
241  event.getByToken(rh_token, cscRecHits);
242 
243  // get CSC segment collection
245  event.getByToken(seg_token, cscSegments);
246 
247  // get the cosmic muons collection
249  if (typeOfSkim == 8) {
250  event.getByToken(sam_token, saMuons);
251  }
252 
253  // get the stand-alone muons collection
256  if (typeOfSkim == 9) {
257  event.getByToken(sam_token, saMuons);
258  event.getByToken(trk_token, tracks);
259  event.getByToken(glm_token, gMuons);
260  }
261 
262  //======================================
263  // evaluate the skimming routines
264  //======================================
265 
266  // basic skimming
267  bool basicEvent = false;
268  if (typeOfSkim == 1 || typeOfSkim == 2) {
269  basicEvent = doCSCSkimming(cscRecHits, cscSegments);
270  }
271 
272  // overlapping chamber skim
273  bool goodOverlapEvent = false;
274  if (typeOfSkim == 3) {
275  goodOverlapEvent = doOverlapSkimming(cscSegments);
276  if (goodOverlapEvent) {
278  }
279  }
280 
281  // messy events skim
282  bool messyEvent = false;
283  if (typeOfSkim == 4) {
284  messyEvent = doMessyEventSkimming(cscRecHits, cscSegments);
285  if (messyEvent) {
286  nEventsMessy++;
287  }
288  }
289 
290  // select events with DIGIs in a certain chamber
291  bool hasChamber = false;
292  if (typeOfSkim == 5) {
293  hasChamber = doCertainChamberSelection(wires, strips);
294  if (hasChamber) {
296  }
297  }
298 
299  // select events in the DT-CSC overlap region
300  bool DTOverlapCandidate = false;
301  if (typeOfSkim == 6) {
302  DTOverlapCandidate = doDTOverlap(cscSegments);
303  if (DTOverlapCandidate) {
305  }
306  }
307 
308  // select halo-like events
309  bool HaloLike = false;
310  if (typeOfSkim == 7) {
311  HaloLike = doHaloLike(cscSegments);
312  if (HaloLike) {
313  nEventsHaloLike++;
314  }
315  }
316 
317  // select long cosmic tracks
318  bool LongSATrack = false;
319  if (typeOfSkim == 8) {
320  LongSATrack = doLongSATrack(saMuons);
321  if (LongSATrack) {
323  }
324  }
325 
326  // select events suitable for a B-field study. They have tracks in the tracker.
327  bool GoodForBFieldStudy = false;
328  if (typeOfSkim == 9) {
329  GoodForBFieldStudy = doBFieldStudySelection(saMuons, tracks, gMuons);
330  if (GoodForBFieldStudy) {
332  }
333  }
334 
335  // set filter flag
336  bool selectThisEvent = false;
337  if (typeOfSkim == 1 || typeOfSkim == 2) {
338  selectThisEvent = basicEvent;
339  }
340  if (typeOfSkim == 3) {
341  selectThisEvent = goodOverlapEvent;
342  }
343  if (typeOfSkim == 4) {
344  selectThisEvent = messyEvent;
345  }
346  if (typeOfSkim == 5) {
347  selectThisEvent = hasChamber;
348  }
349  if (typeOfSkim == 6) {
350  selectThisEvent = DTOverlapCandidate;
351  }
352  if (typeOfSkim == 7) {
353  selectThisEvent = HaloLike;
354  }
355  if (typeOfSkim == 8) {
356  selectThisEvent = LongSATrack;
357  }
358  if (typeOfSkim == 9) {
359  selectThisEvent = GoodForBFieldStudy;
360  }
361 
362  if (selectThisEvent) {
363  nEventsSelected++;
364  }
365 
366  return selectThisEvent;
367 }
368 
369 // ==============================================
370 //
371 // CSC Skimming
372 //
373 // ==============================================
374 
377  // how many RecHits in the collection?
378  int nRecHits = cscRecHits->size();
379 
380  // zero the recHit counter
381  int cntRecHit[600];
382  for (int i = 0; i < 600; i++) {
383  cntRecHit[i] = 0;
384  }
385 
386  // ---------------------
387  // Loop over rechits
388  // ---------------------
389 
391  for (recIt = cscRecHits->begin(); recIt != cscRecHits->end(); recIt++) {
392  // which chamber is it?
393  CSCDetId idrec = (CSCDetId)(*recIt).cscDetId();
394  int kEndcap = idrec.endcap();
395  int kRing = idrec.ring();
396  int kStation = idrec.station();
397  int kChamber = idrec.chamber();
398  int kLayer = idrec.layer();
399 
400  // compute chamber serial number
401  int kSerial = chamberSerial(kEndcap, kStation, kRing, kChamber);
402 
403  // increment recHit counter
404  // (each layer is represented by a different power of 10)
405  int kDigit = (int)pow((float)10., (float)(kLayer - 1));
406  cntRecHit[kSerial] += kDigit;
407 
408  } //end rechit loop
409 
410  // ------------------------------------------------------
411  // Are there chambers with the minimum number of hits?
412  // ------------------------------------------------------
413 
414  int nChambersWithMinimalHits = 0;
415  int nChambersWithMinimalHitsPOS = 0;
416  int nChambersWithMinimalHitsNEG = 0;
417  if (nRecHits > 0) {
418  for (int i = 0; i < 600; i++) {
419  if (cntRecHit[i] > 0) {
420  int nLayersWithHits = 0;
421  float dummy = (float)cntRecHit[i];
422  for (int j = 5; j > -1; j--) {
423  float digit = dummy / pow((float)10., (float)j);
424  int kCount = (int)digit;
425  if (kCount > 0)
426  nLayersWithHits++;
427  dummy = dummy - ((float)kCount) * pow((float)10., (float)j);
428  }
429  if (nLayersWithHits > nLayersWithHitsMinimum) {
430  if (i < 300) {
431  nChambersWithMinimalHitsPOS++;
432  } else {
433  nChambersWithMinimalHitsNEG++;
434  }
435  }
436  }
437  }
438  nChambersWithMinimalHits = nChambersWithMinimalHitsPOS + nChambersWithMinimalHitsNEG;
439  }
440 
441  // how many Segments?
442  int nSegments = cscSegments->size();
443 
444  // ----------------------
445  // fill histograms
446  // ----------------------
447 
448  if (makeHistograms) {
449  hxnRecHits->Fill(nRecHits);
450  if (nRecHits > 0) {
451  hxnSegments->Fill(nSegments);
452  hxnHitChambers->Fill(nChambersWithMinimalHits);
453  }
454  if (nChambersWithMinimalHits > 0) {
455  hxnRecHitsSel->Fill(nRecHits);
456  }
457  }
458 
459  // ----------------------
460  // set the filter flag
461  // ----------------------
462  bool basicEvent = (nChambersWithMinimalHits >= minimumHitChambers) && (nSegments >= minimumSegments);
463 
464  bool chambersOnBothSides =
465  ((nChambersWithMinimalHitsPOS >= minimumHitChambers) && (nChambersWithMinimalHitsNEG >= minimumHitChambers));
466 
467  if (chambersOnBothSides) {
469  }
470 
471  bool selectEvent = false;
472  if (typeOfSkim == 1) {
473  selectEvent = basicEvent;
474  }
475  if (typeOfSkim == 2) {
476  selectEvent = chambersOnBothSides;
477  }
478 
479  // debug
480  LogDebug("[CSCSkim]") << "----- nRecHits = " << nRecHits
481  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits << "\tnSegments = " << nSegments
482  << "\tselect? " << selectEvent << std::endl;
483 
484  /*
485  if ((nChambersWithMinimalHitsPOS >= minimumHitChambers) && (nChambersWithMinimalHitsNEG >= minimumHitChambers)) {
486  std::cout << "\n==========================================================================\n"
487  << "\tinteresting event - chambers hit on both sides\n"
488  << "\t " << nEventsAnalyzed
489  << "\trun " << iRun << "\tevent " << iEvent << std::endl;
490  std::cout << "----- nRecHits = " << nRecHits
491  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits
492  << "\tnSegments = " << nSegments
493  << "\tselect? " << selectEvent << std::endl;
494  for (int i = 0; i < 600; i++) {
495  if (cntRecHit[i] > 0) {
496  cout << "\t\t" << i << "\tcntRecHit= " << cntRecHit[i] << std::endl;
497  }
498  }
499  std::cout << "==========================================================================\n\n" ;
500  }
501  */
502 
503  return selectEvent;
504 }
505 
506 //-------------------------------------------------------------------------------
507 // A module to select events useful in aligning chambers relative to each other
508 // using the overlap regions at the edges (in Xlocal) of the chamber.
509 //-------------------------------------------------------------------------------
511  const int nhitsMinimum = 4;
512  const float chisqMaximum = 100.;
513  const int nAllMaximum = 3;
514 
515  // how many Segments?
516  // int nSegments = cscSegments->size();
517 
518  // zero arrays
519  int nAll[600];
520  int nGood[600];
521  for (int i = 0; i < 600; i++) {
522  nAll[i] = 0;
523  nGood[i] = 0;
524  }
525 
526  // -----------------------
527  // loop over segments
528  // -----------------------
529  for (CSCSegmentCollection::const_iterator it = cscSegments->begin(); it != cscSegments->end(); it++) {
530  // which chamber?
531  CSCDetId id = (CSCDetId)(*it).cscDetId();
532  int kEndcap = id.endcap();
533  int kStation = id.station();
534  int kRing = id.ring();
535  int kChamber = id.chamber();
536  int kSerial = chamberSerial(kEndcap, kStation, kRing, kChamber);
537 
538  // segment information
539  float chisq = (*it).chi2();
540  int nhits = (*it).nRecHits();
541 
542  // is this a good segment?
543  bool goodSegment = (nhits >= nhitsMinimum) && (chisq < chisqMaximum);
544 
545  /*
546  LocalPoint localPos = (*it).localPosition();
547  float segX = localPos.x();
548  float segY = localPos.y();
549  std::cout << "E/S/R/Ch: " << kEndcap << "/" << kStation << "/" << kRing << "/" << kChamber
550  << "\tnhits/chisq: " << nhits << "/" << chisq
551  << "\tX/Y: " << segX << "/" << segY
552  << "\tgood? " << goodSegment << std::endl;
553  */
554 
555  // count
556  nAll[kSerial - 1]++;
557  if (goodSegment)
558  nGood[kSerial]++;
559 
560  } // end loop over segments
561 
562  //----------------------
563  // select the event
564  //----------------------
565 
566  // does any chamber have too many segments?
567  bool messyChamber = false;
568  for (int i = 0; i < 600; i++) {
569  if (nAll[i] > nAllMaximum)
570  messyChamber = true;
571  }
572 
573  // are there consecutive chambers with good segments
574  // (This is a little sloppy but is probably fine for skimming...)
575  bool consecutiveChambers = false;
576  for (int i = 0; i < 599; i++) {
577  if ((nGood[i] > 0) && (nGood[i + 1] > 0))
578  consecutiveChambers = true;
579  }
580 
581  bool selectThisEvent = !messyChamber && consecutiveChambers;
582 
583  return selectThisEvent;
584 }
585 
586 //============================================================
587 //
588 // This module selects events with a large numbere
589 // of recHits and larger number of chambers with hits.
590 //
591 //============================================================
594  // how many RecHits in the collection?
595  int nRecHits = cscRecHits->size();
596 
597  // zero the recHit counter
598  int cntRecHit[600];
599  for (int i = 0; i < 600; i++) {
600  cntRecHit[i] = 0;
601  }
602 
603  // ---------------------
604  // Loop over rechits
605  // ---------------------
606 
608  for (recIt = cscRecHits->begin(); recIt != cscRecHits->end(); recIt++) {
609  // which chamber is it?
610  CSCDetId idrec = (CSCDetId)(*recIt).cscDetId();
611  int kEndcap = idrec.endcap();
612  int kRing = idrec.ring();
613  int kStation = idrec.station();
614  int kChamber = idrec.chamber();
615  int kLayer = idrec.layer();
616 
617  // compute chamber serial number
618  int kSerial = chamberSerial(kEndcap, kStation, kRing, kChamber);
619 
620  // increment recHit counter
621  // (each layer is represented by a different power of 10)
622  int kDigit = (int)pow((float)10., (float)(kLayer - 1));
623  cntRecHit[kSerial] += kDigit;
624 
625  } //end rechit loop
626 
627  // ------------------------------------------------------
628  // Are there chambers with the minimum number of hits?
629  // ------------------------------------------------------
630 
631  int nChambersWithMinimalHits = 0;
632  int nChambersWithMinimalHitsPOS = 0;
633  int nChambersWithMinimalHitsNEG = 0;
634  if (nRecHits > 0) {
635  for (int i = 0; i < 600; i++) {
636  if (cntRecHit[i] > 0) {
637  int nLayersWithHits = 0;
638  float dummy = (float)cntRecHit[i];
639  for (int j = 5; j > -1; j--) {
640  float digit = dummy / pow((float)10., (float)j);
641  int kCount = (int)digit;
642  if (kCount > 0)
643  nLayersWithHits++;
644  dummy = dummy - ((float)kCount) * pow((float)10., (float)j);
645  }
646  if (nLayersWithHits > nLayersWithHitsMinimum) {
647  if (i < 300) {
648  nChambersWithMinimalHitsPOS++;
649  } else {
650  nChambersWithMinimalHitsNEG++;
651  }
652  }
653  }
654  }
655  nChambersWithMinimalHits = nChambersWithMinimalHitsPOS + nChambersWithMinimalHitsNEG;
656  }
657 
658  // how many Segments?
659  int nSegments = cscSegments->size();
660 
661  // ----------------------
662  // fill histograms
663  // ----------------------
664 
666  if (nRecHits > 8) {
667  mevnRecHits0->Fill(nRecHits);
668  mevnChambers0->Fill(nChambersWithMinimalHits);
669  mevnSegments0->Fill(nSegments);
670  }
671  if (nRecHits > 54) {
672  double dummy = (double)nRecHits;
673  if (dummy > 299.9)
674  dummy = 299.9;
675  mevnRecHits1->Fill(dummy);
676  dummy = (double)nChambersWithMinimalHits;
677  if (dummy > 49.9)
678  dummy = 49.9;
679  mevnChambers1->Fill(dummy);
680  dummy = (double)nSegments;
681  if (dummy > 29.9)
682  dummy = 29.9;
683  mevnSegments1->Fill(dummy);
684  }
685  }
686 
687  // ----------------------
688  // set the filter flag
689  // ----------------------
690 
691  bool selectEvent = false;
692  if ((nRecHits > 54) && (nChambersWithMinimalHits > 5)) {
693  selectEvent = true;
694  }
695 
696  // debug
697  LogDebug("[CSCSkim]") << "----- nRecHits = " << nRecHits
698  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits << "\tnSegments = " << nSegments
699  << "\tselect? " << selectEvent << std::endl;
700 
701  /*
702  if (selectEvent) {
703  std::cout << "\n==========================================================================\n"
704  << "\tmessy event!\n"
705  << "\t " << nEventsAnalyzed
706  << "\trun " << iRun << "\tevent " << iEvent << std::endl;
707  std::cout << "----- nRecHits = " << nRecHits
708  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits
709  << "\tnSegments = " << nSegments
710  << "\tselect? " << selectEvent << std::endl;
711  for (int i = 0; i < 600; i++) {
712  if (cntRecHit[i] > 0) {
713  cout << "\t\t" << i << "\tcntRecHit= " << cntRecHit[i] << std::endl;
714  }
715  }
716  std::cout << "==========================================================================\n\n" ;
717  }
718  */
719 
720  return selectEvent;
721 }
722 
723 //============================================================
724 //
725 // Select events with DIGIs are a particular chamber.
726 //
727 //============================================================
730  // Loop through the wire DIGIs, looking for a match
731  bool certainChamberIsPresentInWires = false;
732  for (CSCWireDigiCollection::DigiRangeIterator jw = wires->begin(); jw != wires->end(); jw++) {
733  CSCDetId id = (CSCDetId)(*jw).first;
734  int kEndcap = id.endcap();
735  int kRing = id.ring();
736  int kStation = id.station();
737  int kChamber = id.chamber();
738  if ((kEndcap == whichEndcap) && (kStation == whichStation) && (kRing == whichRing) && (kChamber == whichChamber)) {
739  certainChamberIsPresentInWires = true;
740  }
741  } // end wire loop
742 
743  // Loop through the strip DIGIs, looking for a match
744  bool certainChamberIsPresentInStrips = false;
745  for (CSCStripDigiCollection::DigiRangeIterator js = strips->begin(); js != strips->end(); js++) {
746  CSCDetId id = (CSCDetId)(*js).first;
747  int kEndcap = id.endcap();
748  int kRing = id.ring();
749  int kStation = id.station();
750  int kChamber = id.chamber();
751  if ((kEndcap == whichEndcap) && (kStation == whichStation) && (kRing == whichRing) && (kChamber == whichChamber)) {
752  certainChamberIsPresentInStrips = true;
753  }
754  }
755 
756  bool certainChamberIsPresent = certainChamberIsPresentInWires || certainChamberIsPresentInStrips;
757 
758  return certainChamberIsPresent;
759 }
760 
761 //============================================================
762 //
763 // Select events which *might* probe the DT-CSC overlap region.
764 //
765 //============================================================
767  const float chisqMax = 100.;
768  const int nhitsMin = 5;
769  const int maxNSegments = 3;
770 
771  // initialize
772  bool DTOverlapCandidate = false;
773  int cntMEP13[36];
774  int cntMEN13[36];
775  int cntMEP22[36];
776  int cntMEN22[36];
777  int cntMEP32[36];
778  int cntMEN32[36];
779  for (int i = 0; i < 36; ++i) {
780  cntMEP13[i] = 0;
781  cntMEN13[i] = 0;
782  cntMEP22[i] = 0;
783  cntMEN22[i] = 0;
784  cntMEP32[i] = 0;
785  cntMEN32[i] = 0;
786  }
787 
788  // -----------------------
789  // loop over segments
790  // -----------------------
791 
792  int nSegments = cscSegments->size();
793  if (nSegments < 2)
794  return DTOverlapCandidate;
795 
796  for (CSCSegmentCollection::const_iterator it = cscSegments->begin(); it != cscSegments->end(); it++) {
797  // which chamber?
798  CSCDetId id = (CSCDetId)(*it).cscDetId();
799  int kEndcap = id.endcap();
800  int kStation = id.station();
801  int kRing = id.ring();
802  int kChamber = id.chamber();
803  // segment information
804  float chisq = (*it).chi2();
805  int nhits = (*it).nRecHits();
806  bool goodSegment = (chisq < chisqMax) && (nhits >= nhitsMin);
807  if (goodSegment) {
808  if ((kStation == 1) && (kRing == 3)) {
809  if (kEndcap == 1) {
810  cntMEP13[kChamber - 1]++;
811  }
812  if (kEndcap == 2) {
813  cntMEN13[kChamber - 1]++;
814  }
815  }
816  if ((kStation == 2) && (kRing == 2)) {
817  if (kEndcap == 1) {
818  cntMEP22[kChamber - 1]++;
819  }
820  if (kEndcap == 2) {
821  cntMEN22[kChamber - 1]++;
822  }
823  }
824  if ((kStation == 3) && (kRing == 2)) {
825  if (kEndcap == 1) {
826  cntMEP32[kChamber - 1]++;
827  }
828  if (kEndcap == 2) {
829  cntMEN32[kChamber - 1]++;
830  }
831  }
832  } // this is a good segment
833  } // end loop over segments
834 
835  // ---------------------------------------------
836  // veto messy events
837  // ---------------------------------------------
838  bool tooManySegments = false;
839  for (int i = 0; i < 36; ++i) {
840  if ((cntMEP13[i] > maxNSegments) || (cntMEN13[i] > maxNSegments) || (cntMEP22[i] > maxNSegments) ||
841  (cntMEN22[i] > maxNSegments) || (cntMEP32[i] > maxNSegments) || (cntMEN32[i] > maxNSegments))
842  tooManySegments = true;
843  }
844  if (tooManySegments) {
845  return DTOverlapCandidate;
846  }
847 
848  // ---------------------------------------------
849  // check for relevant matchup of segments
850  // ---------------------------------------------
851  bool matchup = false;
852  for (int i = 0; i < 36; ++i) {
853  if ((cntMEP13[i] > 0) && (cntMEP22[i] + cntMEP32[i] > 0)) {
854  matchup = true;
855  }
856  if ((cntMEN13[i] > 0) && (cntMEN22[i] + cntMEN32[i] > 0)) {
857  matchup = true;
858  }
859  }
860  /*
861  if (matchup) {
862  std::cout << "\tYYY looks like a good event. Select!\n";
863  std::cout << "-- pos endcap --\n"
864  << "ME1/3: ";
865  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP13[k];}
866  std::cout << "\nME2/2: ";
867  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP22[k];}
868  std::cout << "\nME3/2: ";
869  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP32[k];}
870  std::cout << std::endl;
871  }
872  */
873 
874  // set the selection flag
875  DTOverlapCandidate = matchup;
876  return DTOverlapCandidate;
877 }
878 
879 //============================================================
880 //
881 // Select events which register in the inner parts of
882 // stations 1, 2, 3 and 4.
883 //
884 //============================================================
886  const float chisqMax = 100.;
887  const int nhitsMin = 5; // on a segment
888  const int maxNSegments = 3; // in a chamber
889 
890  // initialize
891  bool HaloLike = false;
892  int cntMEP11[36];
893  int cntMEN11[36];
894  int cntMEP12[36];
895  int cntMEN12[36];
896  int cntMEP21[36];
897  int cntMEN21[36];
898  int cntMEP31[36];
899  int cntMEN31[36];
900  int cntMEP41[36];
901  int cntMEN41[36];
902  for (int i = 0; i < 36; ++i) {
903  cntMEP11[i] = 0;
904  cntMEN11[i] = 0;
905  cntMEP12[i] = 0;
906  cntMEN12[i] = 0;
907  cntMEP21[i] = 0;
908  cntMEN21[i] = 0;
909  cntMEP31[i] = 0;
910  cntMEN31[i] = 0;
911  cntMEP41[i] = 0;
912  cntMEN41[i] = 0;
913  }
914 
915  // -----------------------
916  // loop over segments
917  // -----------------------
918  int nSegments = cscSegments->size();
919  if (nSegments < 4)
920  return HaloLike;
921 
922  for (CSCSegmentCollection::const_iterator it = cscSegments->begin(); it != cscSegments->end(); it++) {
923  // which chamber?
924  CSCDetId id = (CSCDetId)(*it).cscDetId();
925  int kEndcap = id.endcap();
926  int kStation = id.station();
927  int kRing = id.ring();
928  int kChamber = id.chamber();
929  // segment information
930  float chisq = (*it).chi2();
931  int nhits = (*it).nRecHits();
932  bool goodSegment = (chisq < chisqMax) && (nhits >= nhitsMin);
933  if (goodSegment) {
934  if ((kStation == 1) && (kRing == 1)) {
935  if (kEndcap == 1) {
936  cntMEP11[kChamber - 1]++;
937  }
938  if (kEndcap == 2) {
939  cntMEN11[kChamber - 1]++;
940  }
941  }
942  if ((kStation == 1) && (kRing == 2)) {
943  if (kEndcap == 1) {
944  cntMEP12[kChamber - 1]++;
945  }
946  if (kEndcap == 2) {
947  cntMEN12[kChamber - 1]++;
948  }
949  }
950  if ((kStation == 2) && (kRing == 1)) {
951  if (kEndcap == 1) {
952  cntMEP21[kChamber - 1]++;
953  }
954  if (kEndcap == 2) {
955  cntMEN21[kChamber - 1]++;
956  }
957  }
958  if ((kStation == 3) && (kRing == 1)) {
959  if (kEndcap == 1) {
960  cntMEP31[kChamber - 1]++;
961  }
962  if (kEndcap == 2) {
963  cntMEN31[kChamber - 1]++;
964  }
965  }
966  if ((kStation == 4) && (kRing == 1)) {
967  if (kEndcap == 1) {
968  cntMEP41[kChamber - 1]++;
969  }
970  if (kEndcap == 2) {
971  cntMEN41[kChamber - 1]++;
972  }
973  }
974  } // this is a good segment
975  } // end loop over segments
976 
977  // ---------------------------------------------
978  // veto messy events
979  // ---------------------------------------------
980  bool tooManySegments = false;
981  for (int i = 0; i < 36; ++i) {
982  if ((cntMEP11[i] > 3 * maxNSegments) || (cntMEN11[i] > 3 * maxNSegments) || (cntMEP12[i] > maxNSegments) ||
983  (cntMEN12[i] > maxNSegments) || (cntMEP21[i] > maxNSegments) || (cntMEN21[i] > maxNSegments) ||
984  (cntMEP31[i] > maxNSegments) || (cntMEN31[i] > maxNSegments) || (cntMEP41[i] > maxNSegments) ||
985  (cntMEN41[i] > maxNSegments))
986  tooManySegments = true;
987  }
988  if (tooManySegments) {
989  return HaloLike;
990  }
991 
992  // ---------------------------------------------
993  // check for relevant matchup of segments
994  // ---------------------------------------------
995  bool matchup = false;
996  for (int i = 0; i < 36; ++i) {
997  if ((cntMEP11[i] + cntMEP12[i] > 0) && (cntMEP21[i] > 0) && (cntMEP31[i] > 0) && (cntMEP41[i] > 0)) {
998  matchup = true;
999  }
1000  if ((cntMEN11[i] + cntMEN12[i] > 0) && (cntMEN21[i] > 0) && (cntMEN31[i] > 0) && (cntMEN41[i] > 0)) {
1001  matchup = true;
1002  }
1003  }
1004  /*
1005  if (matchup) {
1006  std::cout << "\tYYY looks like a good event. Select!\n";
1007  std::cout << "-- pos endcap --\n"
1008  << "ME1/1: ";
1009  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP11[k];}
1010  std::cout << "\nME1/2: ";
1011  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP12[k];}
1012  std::cout << "\nME2/1: ";
1013  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP21[k];}
1014  std::cout << "\nME3/1: ";
1015  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP31[k];}
1016  std::cout << "\nME4/1: ";
1017  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP41[k];}
1018  std::cout << std::endl;
1019  std::cout << "-- neg endcap --\n"
1020  << "ME1/1: ";
1021  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN11[k];}
1022  std::cout << "\nME1/2: ";
1023  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN12[k];}
1024  std::cout << "\nME2/1: ";
1025  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN21[k];}
1026  std::cout << "\nME3/1: ";
1027  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN31[k];}
1028  std::cout << "\nME4/1: ";
1029  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN41[k];}
1030  std::cout << std::endl;
1031  std::cout << "\tn Analyzed = " << nEventsAnalyzed << "\tn Halo-like = " << nEventsHaloLike << std::endl;
1032  }
1033  */
1034 
1035  // set the selection flag
1036  HaloLike = matchup;
1037  return HaloLike;
1038 }
1039 
1040 //--------------------------------------------------------------
1041 // select events with at least one "long" stand-alone muon
1042 //--------------------------------------------------------------
1044  const float zDistanceMax = 2500.;
1045  const float zDistanceMin = 700.;
1046  const int nCSCHitsMin = 25;
1047  const int nCSCHitsMax = 50;
1048  const float zInnerMax = 80000.;
1049 
1050  const int nNiceMuonsMin = 1;
1051 
1052  //
1053  // Loop through the track collection and test each one
1054  //
1055 
1056  int nNiceMuons = 0;
1057 
1058  for (reco::TrackCollection::const_iterator muon = saMuons->begin(); muon != saMuons->end(); ++muon) {
1059  // basic information
1060  math::XYZVector innerMo = muon->innerMomentum();
1061  GlobalVector im(innerMo.x(), innerMo.y(), innerMo.z());
1062  math::XYZPoint innerPo = muon->innerPosition();
1063  GlobalPoint ip(innerPo.x(), innerPo.y(), innerPo.z());
1064  math::XYZPoint outerPo = muon->outerPosition();
1065  GlobalPoint op(outerPo.x(), outerPo.y(), outerPo.z());
1066  float zInner = ip.z();
1067  float zOuter = op.z();
1068  float zDistance = fabs(zOuter - zInner);
1069 
1070  // loop over hits
1071  int nDTHits = 0;
1072  int nCSCHits = 0;
1073  for (trackingRecHit_iterator hit = muon->recHitsBegin(); hit != muon->recHitsEnd(); ++hit) {
1074  const DetId detId((*hit)->geographicalId());
1075  if (detId.det() == DetId::Muon) {
1076  if (detId.subdetId() == MuonSubdetId::DT) {
1077  //DTChamberId dtId(detId.rawId());
1078  //int chamberId = dtId.sector();
1079  nDTHits++;
1080  } else if (detId.subdetId() == MuonSubdetId::CSC) {
1081  //CSCDetId cscId(detId.rawId());
1082  //int chamberId = cscId.chamber();
1083  nCSCHits++;
1084  }
1085  }
1086  }
1087 
1088  // is this a nice muon?
1089  if ((zDistance < zDistanceMax) && (zDistance > zDistanceMin) && (nCSCHits > nCSCHitsMin) &&
1090  (nCSCHits < nCSCHitsMax) && (min(fabs(zInner), fabs(zOuter)) < zInnerMax) &&
1091  (fabs(innerMo.z()) > 0.000000001)) {
1092  nNiceMuons++;
1093  }
1094  }
1095 
1096  bool select = (nNiceMuons >= nNiceMuonsMin);
1097 
1098  return select;
1099 }
1100 
1101 //============================================================
1102 //
1103 // Select events which are good for B-field studies.
1104 //
1105 // These events have a good track in the tracker.
1106 //
1107 // D.Dibur and M.Schmitt
1108 //============================================================
1112  bool acceptThisEvent = false;
1113 
1114  //-----------------------------------
1115  // examine the stand-alone tracks
1116  //-----------------------------------
1117  int nGoodSAMuons = 0;
1118  for (reco::TrackCollection::const_iterator muon = saMuons->begin(); muon != saMuons->end(); ++muon) {
1119  float preco = muon->p();
1120 
1121  math::XYZPoint innerPo = muon->innerPosition();
1122  GlobalPoint iPnt(innerPo.x(), innerPo.y(), innerPo.z());
1123  math::XYZPoint outerPo = muon->outerPosition();
1124  GlobalPoint oPnt(outerPo.x(), outerPo.y(), outerPo.z());
1125  float zLength = abs(iPnt.z() - oPnt.z());
1126 
1127  math::XYZVector innerMom = muon->innerMomentum();
1128  GlobalVector iP(innerMom.x(), innerMom.y(), innerMom.z());
1129  math::XYZVector outerMom = muon->outerMomentum();
1130  GlobalVector oP(outerMom.x(), outerMom.y(), outerMom.z());
1131 
1132  const float zRef = 300.;
1133  float xExt = 10000.;
1134  float yExt = 10000.;
1135  if (abs(oPnt.z()) < abs(iPnt.z())) {
1136  float deltaZ = 0.;
1137  if (oPnt.z() > 0) {
1138  deltaZ = zRef - oPnt.z();
1139  } else {
1140  deltaZ = -zRef - oPnt.z();
1141  }
1142  xExt = oPnt.x() + deltaZ * oP.x() / oP.z();
1143  yExt = oPnt.y() + deltaZ * oP.y() / oP.z();
1144  } else {
1145  float deltaZ = 0.;
1146  if (iPnt.z() > 0) {
1147  deltaZ = zRef - iPnt.z();
1148  } else {
1149  deltaZ = -zRef - iPnt.z();
1150  }
1151  xExt = iPnt.x() + deltaZ * iP.x() / iP.z();
1152  yExt = iPnt.y() + deltaZ * iP.y() / iP.z();
1153  }
1154  float rExt = sqrt(xExt * xExt + yExt * yExt);
1155 
1156  int kHit = 0;
1157  int nDTHits = 0;
1158  int nCSCHits = 0;
1159  for (trackingRecHit_iterator hit = muon->recHitsBegin(); hit != muon->recHitsEnd(); ++hit) {
1160  ++kHit;
1161  const DetId detId((*hit)->geographicalId());
1162  if (detId.det() == DetId::Muon) {
1163  if (detId.subdetId() == MuonSubdetId::DT) {
1164  nDTHits++;
1165  } else if (detId.subdetId() == MuonSubdetId::CSC) {
1166  nCSCHits++;
1167  }
1168  }
1169  } // end loop over hits
1170 
1171  float zInner = -1.;
1172  if (nCSCHits >= nCSCHitsMin) {
1173  if (abs(iPnt.z()) < abs(iPnt.z())) {
1174  zInner = iPnt.z();
1175  } else {
1176  zInner = oPnt.z();
1177  }
1178  }
1179 
1180  bool goodSAMuon = (preco > pMin) && (zLength > zLengthMin) && (nCSCHits >= nCSCHitsMin) && (zInner < zInnerMax) &&
1181  (rExt < rExtMax);
1182 
1183  if (goodSAMuon) {
1184  nGoodSAMuons++;
1185  }
1186 
1187  } // end loop over stand-alone muon collection
1188 
1189  //-----------------------------------
1190  // examine the tracker tracks
1191  //-----------------------------------
1192  int nGoodTracks = 0;
1193  for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
1194  float preco = track->p();
1195  int n = track->recHitsSize();
1196 
1197  math::XYZPoint innerPo = track->innerPosition();
1198  GlobalPoint iPnt(innerPo.x(), innerPo.y(), innerPo.z());
1199  math::XYZPoint outerPo = track->outerPosition();
1200  GlobalPoint oPnt(outerPo.x(), outerPo.y(), outerPo.z());
1201  float zLength = abs(iPnt.z() - oPnt.z());
1202 
1203  math::XYZVector innerMom = track->innerMomentum();
1204  GlobalVector iP(innerMom.x(), innerMom.y(), innerMom.z());
1205  math::XYZVector outerMom = track->outerMomentum();
1206  GlobalVector oP(outerMom.x(), outerMom.y(), outerMom.z());
1207 
1208  const float zRef = 300.;
1209  float xExt = 10000.;
1210  float yExt = 10000.;
1211  if (abs(oPnt.z()) > abs(iPnt.z())) {
1212  float deltaZ = 0.;
1213  if (oPnt.z() > 0) {
1214  deltaZ = zRef - oPnt.z();
1215  } else {
1216  deltaZ = -zRef - oPnt.z();
1217  }
1218  xExt = oPnt.x() + deltaZ * oP.x() / oP.z();
1219  yExt = oPnt.y() + deltaZ * oP.y() / oP.z();
1220  } else {
1221  float deltaZ = 0.;
1222  if (iPnt.z() > 0) {
1223  deltaZ = zRef - iPnt.z();
1224  } else {
1225  deltaZ = -zRef - iPnt.z();
1226  }
1227  xExt = iPnt.x() + deltaZ * iP.x() / iP.z();
1228  yExt = iPnt.y() + deltaZ * iP.y() / iP.z();
1229  }
1230  float rExt = sqrt(xExt * xExt + yExt * yExt);
1231 
1232  bool goodTrack = (preco > pMin) && (n >= nTrHitsMin) && (zLength > zLengthTrMin) && (rExt < rExtMax);
1233 
1234  if (goodTrack) {
1235  nGoodTracks++;
1236  }
1237 
1238  } // end loop over tracker tracks
1239 
1240  //-----------------------------------
1241  // examine the global muons
1242  //-----------------------------------
1243  int nGoodGlobalMuons = 0;
1244  for (reco::MuonCollection::const_iterator global = gMuons->begin(); global != gMuons->end(); ++global) {
1245  if (global->isGlobalMuon()) {
1246  float pDef = global->p();
1247  float redChiSq = global->globalTrack()->normalizedChi2();
1248  const reco::HitPattern& hp = (global->globalTrack())->hitPattern();
1249  // int nTotalHits = hp.numberOfHits();
1250  // int nValidHits = hp.numberOfValidHits();
1251  int nTrackerHits = hp.numberOfValidTrackerHits();
1252  // int nPixelHits = hp.numberOfValidPixelHits();
1253  // int nStripHits = hp.numberOfValidStripHits();
1254 
1255  int nDTHits = 0;
1256  int nCSCHits = 0;
1257  for (trackingRecHit_iterator hit = (global->globalTrack())->recHitsBegin();
1258  hit != (global->globalTrack())->recHitsEnd();
1259  ++hit) {
1260  const DetId detId((*hit)->geographicalId());
1261  if (detId.det() == DetId::Muon) {
1262  if (detId.subdetId() == MuonSubdetId::DT) {
1263  nDTHits++;
1264  } else if (detId.subdetId() == MuonSubdetId::CSC) {
1265  nCSCHits++;
1266  }
1267  }
1268  } // end loop over hits
1269 
1270  bool goodGlobalMuon =
1271  (pDef > pMin) && (nTrackerHits >= nValidHitsMin) && (nCSCHits >= nCSCHitsMin) && (redChiSq < redChiSqMax);
1272 
1273  if (goodGlobalMuon) {
1274  nGoodGlobalMuons++;
1275  }
1276 
1277  } // this is a global muon
1278  } // end loop over stand-alone muon collection
1279 
1280  //-----------------------------------
1281  // do we accept this event?
1282  //-----------------------------------
1283 
1284  acceptThisEvent = ((nGoodSAMuons > 0) && (nGoodTracks > 0)) || (nGoodGlobalMuons > 0);
1285 
1286  return acceptThisEvent;
1287 }
1288 
1289 //--------------------------------------------------------------
1290 // Compute a serial number for the chamber.
1291 // This is useful when filling histograms and working with arrays.
1292 //--------------------------------------------------------------
1293 int CSCSkim::chamberSerial(int kEndcap, int kStation, int kRing, int kChamber) {
1294  int kSerial = kChamber;
1295  if (kStation == 1 && kRing == 1) {
1296  kSerial = kChamber;
1297  }
1298  if (kStation == 1 && kRing == 2) {
1299  kSerial = kChamber + 36;
1300  }
1301  if (kStation == 1 && kRing == 3) {
1302  kSerial = kChamber + 72;
1303  }
1304  if (kStation == 1 && kRing == 4) {
1305  kSerial = kChamber;
1306  }
1307  if (kStation == 2 && kRing == 1) {
1308  kSerial = kChamber + 108;
1309  }
1310  if (kStation == 2 && kRing == 2) {
1311  kSerial = kChamber + 126;
1312  }
1313  if (kStation == 3 && kRing == 1) {
1314  kSerial = kChamber + 162;
1315  }
1316  if (kStation == 3 && kRing == 2) {
1317  kSerial = kChamber + 180;
1318  }
1319  if (kStation == 4 && kRing == 1) {
1320  kSerial = kChamber + 216;
1321  }
1322  if (kStation == 4 && kRing == 2) {
1323  kSerial = kChamber + 234;
1324  } // one day...
1325  if (kEndcap == 2) {
1326  kSerial = kSerial + 300;
1327  }
1328  return kSerial;
1329 }
1330 
1331 //define this as a plug-in
edm::EDGetTokenT< CSCSegmentCollection > seg_token
Definition: CSCSkim.h:149
TH1F * hxnRecHitsSel
Definition: CSCSkim.h:182
edm::EDGetTokenT< reco::TrackCollection > sam_token
Definition: CSCSkim.h:150
float zInnerMax
Definition: CSCSkim.h:171
float rExtMax
Definition: CSCSkim.h:174
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
bool makeHistogramsForMessyEvents
Definition: CSCSkim.h:162
int nEventsForBFieldStudies
Definition: CSCSkim.h:126
int whichEndcap
Definition: CSCSkim.h:163
int nCSCHitsMin
Definition: CSCSkim.h:170
int nLayersWithHitsMinimum
Definition: CSCSkim.h:157
bool filter(edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: CSCSkim.cc:215
float zLengthMin
Definition: CSCSkim.h:169
int nEventsChambersBothSides
Definition: CSCSkim.h:119
TH1F * mevnRecHits0
Definition: CSCSkim.h:183
TH1F * mevnRecHits1
Definition: CSCSkim.h:186
T z() const
Definition: PV3DBase.h:61
TH1F * mevnSegments0
Definition: CSCSkim.h:185
~CSCSkim() override
Definition: CSCSkim.cc:114
edm::EDGetTokenT< CSCWireDigiCollection > wds_token
Definition: CSCSkim.h:143
bool doLongSATrack(edm::Handle< reco::TrackCollection > saTracks)
Definition: CSCSkim.cc:1043
TH1F * mevnSegments1
Definition: CSCSkim.h:188
std::string outputFileName
Definition: CSCSkim.h:136
TH1F * hxnRecHits
Definition: CSCSkim.h:179
float pMin
Definition: CSCSkim.h:168
TH1F * xxnCSCHits
Definition: CSCSkim.h:190
int layer() const
Definition: CSCDetId.h:56
int minimumHitChambers
Definition: CSCSkim.h:158
bool doOverlapSkimming(edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:510
int nEventsDTOverlap
Definition: CSCSkim.h:123
bool doDTOverlap(edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:766
TH1F * xxnTrackerHits
Definition: CSCSkim.h:190
int nEventsOverlappingChambers
Definition: CSCSkim.h:120
int typeOfSkim
Definition: CSCSkim.h:156
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
edm::EDGetTokenT< reco::TrackCollection > trk_token
Definition: CSCSkim.h:151
int nEventsSelected
Definition: CSCSkim.h:118
float zLengthTrMin
Definition: CSCSkim.h:173
int whichChamber
Definition: CSCSkim.h:166
bool goodTrack(const reco::Track *pTrack, math::XYZPoint leadPV, trackSelectionParameters parameters, bool debug=false)
TH1F * hxnSegments
Definition: CSCSkim.h:180
int minimumSegments
Definition: CSCSkim.h:159
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< reco::MuonCollection > glm_token
Definition: CSCSkim.h:152
std::string histogramFileName
Definition: CSCSkim.h:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< CSCStripDigiCollection > sds_token
Definition: CSCSkim.h:144
int nValidHitsMin
Definition: CSCSkim.h:176
int chamber() const
Definition: CSCDetId.h:62
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< CSCRecHit2DCollection > rh_token
Definition: CSCSkim.h:148
bool makeHistograms
Definition: CSCSkim.h:161
int nEventsLongSATrack
Definition: CSCSkim.h:125
static const std::string kLayer("layer")
int nTrHitsMin
Definition: CSCSkim.h:172
bool doCSCSkimming(edm::Handle< CSCRecHit2DCollection > cscRecHits, edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:375
bool doCertainChamberSelection(edm::Handle< CSCWireDigiCollection > wires, edm::Handle< CSCStripDigiCollection > strips)
Definition: CSCSkim.cc:728
bool doBFieldStudySelection(edm::Handle< reco::TrackCollection > saTracks, edm::Handle< reco::TrackCollection > Tracks, edm::Handle< reco::MuonCollection > gMuons)
Definition: CSCSkim.cc:1109
int nEventsAnalyzed
Definition: CSCSkim.h:117
Log< level::Info, false > LogInfo
int nEventsCertainChamber
Definition: CSCSkim.h:122
Definition: DetId.h:17
auto const & tracks
cannot be loose
int station() const
Definition: CSCDetId.h:79
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void beginJob() override
Definition: CSCSkim.cc:119
int endcap() const
Definition: CSCDetId.h:85
TH1F * mevnChambers1
Definition: CSCSkim.h:187
TFile * theHistogramFile
Definition: CSCSkim.h:133
int whichStation
Definition: CSCSkim.h:164
bool doMessyEventSkimming(edm::Handle< CSCRecHit2DCollection > cscRecHits, edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:592
HLT enums.
int nEventsMessy
Definition: CSCSkim.h:121
int whichRing
Definition: CSCSkim.h:165
const edm::ESGetToken< CSCGeometry, MuonGeometryRecord > m_CSCGeomToken
Definition: CSCSkim.h:140
TH1F * hxnHitChambers
Definition: CSCSkim.h:181
strips
#turn off noise in all subdetectors simHcalUnsuppressedDigis.doNoise = False mix.digitizers.hcal.doNoise = False simEcalUnsuppressedDigis.doNoise = False mix.digitizers.ecal.doNoise = False simEcalUnsuppressedDigis.doESNoise = False simSiPixelDigis.AddNoise = False mix.digitizers.pixel.AddNoise = False simSiStripDigis.Noise = False mix.digitizers.strip.AddNoise = False
Definition: DigiDM_cff.py:32
int chamberSerial(int kE, int kS, int kR, int kCh)
Definition: CSCSkim.cc:1293
int iEvent
Definition: CSCSkim.h:130
static constexpr int DT
Definition: MuonSubdetId.h:11
int ring() const
Definition: CSCDetId.h:68
void endJob() override
Definition: CSCSkim.cc:167
CSCSkim(const edm::ParameterSet &pset)
Definition: CSCSkim.cc:51
edm::EDGetTokenT< CSCWireDigiCollection > wdr_token
Definition: CSCSkim.h:145
int nEventsHaloLike
Definition: CSCSkim.h:124
static constexpr int CSC
Definition: MuonSubdetId.h:12
edm::EDGetTokenT< CSCStripDigiCollection > sdr_token
Definition: CSCSkim.h:146
int iRun
Definition: CSCSkim.h:129
TH1F * xxnValidHits
Definition: CSCSkim.h:190
TH1F * mevnChambers0
Definition: CSCSkim.h:184
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
TH1F * xxredChiSq
Definition: CSCSkim.h:190
Definition: event.py:1
bool demandChambersBothSides
Definition: CSCSkim.h:160
#define LogDebug(id)
TH1F * xxP
Definition: CSCSkim.h:190
bool doHaloLike(edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:885
float redChiSqMax
Definition: CSCSkim.h:175