CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/DPGAnalysis/Skims/src/CSCSkim.cc

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