CMS 3D CMS Logo

CSCDigiSuppressor.cc

Go to the documentation of this file.
00001 #include "L1Trigger/CSCTriggerPrimitives/src/CSCDigiSuppressor.h"
00002 #include "DataFormats/Common/interface/Handle.h"
00003 #include <DataFormats/CSCDigi/interface/CSCCorrelatedLCTDigiCollection.h>
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/Framework/interface/EventSetup.h"
00006 
00007 
00008 CSCDigiSuppressor::CSCDigiSuppressor(const edm::ParameterSet& ps)
00009 : theLCTTag(ps.getParameter<edm::InputTag>("lctTag")),
00010   theStripDigiTag(ps.getParameter<edm::InputTag>("stripDigiTag"))
00011 {
00012   produces<CSCStripDigiCollection>("MuonCSCSuppressedStripDigi");
00013 }
00014 
00015 
00016 void CSCDigiSuppressor::produce(edm::Event& e, const edm::EventSetup& eventSetup) 
00017 {
00018   edm::Handle<CSCStripDigiCollection> oldStripDigis;
00019   e.getByLabel(theStripDigiTag, oldStripDigis);
00020   if (!oldStripDigis.isValid()) {
00021     edm::LogError("CSCDigiValidation") << "Cannot get strips by label "
00022                                        << theStripDigiTag.encode();
00023   }
00024 
00025   edm::Handle<CSCCorrelatedLCTDigiCollection> lcts;
00026   e.getByLabel(theLCTTag, lcts);
00027 
00028   std::auto_ptr<CSCStripDigiCollection> newStripDigis(new CSCStripDigiCollection());
00029 
00030 
00031   for(CSCCorrelatedLCTDigiCollection::DigiRangeIterator detUnitItr = lcts->begin();
00032       detUnitItr != lcts->end(); ++detUnitItr)
00033   {
00034     const CSCDetId& id = (*detUnitItr).first;
00035     const CSCCorrelatedLCTDigiCollection::Range& range = (*detUnitItr).second;
00036     std::list<int> keyStrips;
00037     for (CSCCorrelatedLCTDigiCollection::const_iterator digiItr = range.first;
00038          digiItr != range.second; digiItr++) 
00039     {
00040       // convert from 0..159 to 1..80
00041       keyStrips.push_back(digiItr->getStrip()/2+1);
00042     }
00043 
00044     fillDigis(id, keyStrips, *oldStripDigis, *newStripDigis);
00045   }
00046 
00047   e.put(newStripDigis, "MuonCSCSuppressedStripDigi");
00048 }
00049 
00050 
00051 void CSCDigiSuppressor::fillDigis(const CSCDetId & id, const std::list<int> & keyStrips, 
00052                                        const CSCStripDigiCollection & oldStripDigis,
00053                                        CSCStripDigiCollection & newStripDigis)
00054 {
00055   std::list<int> cfebs = cfebsToRead(id, keyStrips);
00056   CSCStripDigiCollection::Range chamberDigis = oldStripDigis.get(id);
00057   // strips are sorted by layer
00058   for(int layer = 1; layer <= 6; ++layer)
00059   {
00060     CSCDetId layerId(id.rawId()+layer);
00061     for(CSCStripDigiCollection::const_iterator digiItr = chamberDigis.first;
00062         digiItr != chamberDigis.second; ++digiItr)
00063     {
00064       int cfeb = (digiItr->getStrip()-1)/16; 
00065       if(std::find(cfebs.begin(), cfebs.end(), cfeb) != cfebs.end())
00066       {
00067         newStripDigis.insertDigi(layerId, *digiItr);
00068       }
00069     }
00070   }
00071 }
00072 
00073 
00074 std::list<int>
00075 CSCDigiSuppressor::cfebsToRead(const CSCDetId & id, const std::list<int> & keyStrips) const
00076 {
00077   // always accept ME1A, because it's too much trouble looking
00078   // for LCTs in ME11
00079   if(id.station() == 1 && id.ring() == 4)
00080   {
00081     return std::list<int>(1, 0.);
00082   }
00083 
00084   int maxCFEBs = (id.station() == 1) ? 4 : 5;
00085   if(id.station() == 1 && id.ring() == 2) maxCFEBs = 5;
00086 
00087   //copied from CSCStripElectronicsSim
00088   std::list<int> cfebs;
00089   for(std::list<int>::const_iterator keyStripItr = keyStrips.begin(); 
00090       keyStripItr != keyStrips.end(); ++keyStripItr)
00091   {
00092     int cfeb = ((*keyStripItr)-1)/16;
00093     cfebs.push_back(cfeb);
00094     int remainder = ((*keyStripItr)-1)%16;
00095     // if we're within 3 strips of an edge, take neighboring CFEB, too
00096     if(remainder <= 2 && cfeb != 0)
00097     {
00098       cfebs.push_back(cfeb-1);
00099     }
00100 
00101     if(remainder >= 13 && cfeb < maxCFEBs) 
00102     {
00103       cfebs.push_back(cfeb+1);
00104     }
00105   }
00106   cfebs.sort();
00107   cfebs.unique();
00108   return cfebs;
00109 }
00110 
00111 
00112 

Generated on Tue Jun 9 17:39:59 2009 for CMSSW by  doxygen 1.5.4