Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:43

0001 // -*- C++ -*-
0002 //
0003 // Class:      HLTRPCFilter
0004 //
0005 /**\class RPCPathChambFilter HLTRPCFilter.cc HLTRPCFilter.cc
0006 
0007  Description: <one line class summary>
0008 
0009  Implementation:
0010      <Notes on implementation>
0011 */
0012 //
0013 // Original Author:  Camilo Andres Carrillo Montoya
0014 //         Created:  Thu Oct 29 11:04:22 CET 2009
0015 //
0016 //
0017 
0018 #include "HLTRPCFilter.h"
0019 
0020 // system include files
0021 #include <memory>
0022 
0023 //
0024 // constants, enums and typedefs
0025 //
0026 
0027 //
0028 // static data member definitions
0029 //
0030 
0031 //
0032 // constructors and destructor
0033 //
0034 
0035 HLTRPCFilter::HLTRPCFilter(const edm::ParameterSet& config)
0036     : rpcRecHitsToken(consumes<RPCRecHitCollection>(config.getParameter<edm::InputTag>("rpcRecHits"))),
0037       rpcDTPointsToken(consumes<RPCRecHitCollection>(config.getParameter<edm::InputTag>("rpcDTPoints"))),
0038       rpcCSCPointsToken(consumes<RPCRecHitCollection>(config.getParameter<edm::InputTag>("rpcCSCPoints"))),
0039       rangestrips(config.getUntrackedParameter<double>("rangestrips", 1.)) {}
0040 
0041 HLTRPCFilter::~HLTRPCFilter() {
0042   // do anything here that needs to be done at desctruction time
0043   // (e.g. close files, deallocate resources etc.)
0044 }
0045 
0046 void HLTRPCFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0047   edm::ParameterSetDescription desc;
0048   desc.add<edm::InputTag>("rpcRecHits", edm::InputTag("hltRpcRecHits"));
0049   desc.add<edm::InputTag>("rpcDTPoints", edm::InputTag("rpcPointProducer", "RPCDTExtrapolatedPoints"));
0050   desc.add<edm::InputTag>("rpcCSCPoints", edm::InputTag("rpcPointProducer", "RPCCSCExtrapolatedPoints"));
0051   desc.addUntracked<double>("rangestrips", 4.0);
0052   descriptions.add("hltRPCFilter", desc);
0053 }
0054 
0055 //
0056 // member functions
0057 //
0058 
0059 // ------------ method called on each new Event  ------------
0060 bool HLTRPCFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0061   edm::Handle<RPCRecHitCollection> rpcHits;
0062   iEvent.getByToken(rpcRecHitsToken, rpcHits);
0063 
0064   RPCRecHitCollection::const_iterator rpcPoint;
0065 
0066   if (rpcHits->begin() == rpcHits->end()) {
0067     //std::cout<<" skipped preventing no RPC runs"<<std::endl;
0068     return false;
0069   }
0070 
0071   edm::Handle<RPCRecHitCollection> rpcDTPoints;
0072   iEvent.getByToken(rpcDTPointsToken, rpcDTPoints);
0073 
0074   edm::Handle<RPCRecHitCollection> rpcCSCPoints;
0075   iEvent.getByToken(rpcCSCPointsToken, rpcCSCPoints);
0076 
0077   float cluSize = 0;
0078 
0079   //DTPart
0080 
0081   for (rpcPoint = rpcDTPoints->begin(); rpcPoint != rpcDTPoints->end(); rpcPoint++) {
0082     LocalPoint PointExtrapolatedRPCFrame = rpcPoint->localPosition();
0083     RPCDetId rpcId = rpcPoint->rpcId();
0084     typedef std::pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
0085     rangeRecHits recHitCollection = rpcHits->get(rpcId);
0086     if (recHitCollection.first == recHitCollection.second) {
0087       //std::cout<<"DT passed, no rechits for this RPCId =  "<<rpcId<<std::endl;
0088       return true;
0089     }
0090     float minres = 3000.;
0091     RPCRecHitCollection::const_iterator recHit;
0092     for (recHit = recHitCollection.first; recHit != recHitCollection.second; recHit++) {
0093       LocalPoint recHitPos = recHit->localPosition();
0094       float res = PointExtrapolatedRPCFrame.x() - recHitPos.x();
0095       if (fabs(res) < fabs(minres)) {
0096         minres = res;
0097         cluSize = recHit->clusterSize();
0098       }
0099     }
0100     if (fabs(minres) >= (rangestrips + cluSize * 0.5) * 3) {  //3 is a typyical strip width for RPCs
0101       //std::cout<<"DT passed, RecHits but far away "<<rpcId<<std::endl;
0102       return true;
0103     }
0104   }
0105 
0106   //CSCPart
0107 
0108   for (rpcPoint = rpcCSCPoints->begin(); rpcPoint != rpcCSCPoints->end(); rpcPoint++) {
0109     LocalPoint PointExtrapolatedRPCFrame = rpcPoint->localPosition();
0110     RPCDetId rpcId = rpcPoint->rpcId();
0111     typedef std::pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
0112     rangeRecHits recHitCollection = rpcHits->get(rpcId);
0113     if (recHitCollection.first == recHitCollection.second) {
0114       //std::cout<<"CSC passed, no rechits for this RPCId =  "<<rpcId<<std::endl;
0115       return true;
0116     }
0117     float minres = 3000.;
0118     RPCRecHitCollection::const_iterator recHit;
0119     for (recHit = recHitCollection.first; recHit != recHitCollection.second; recHit++) {
0120       LocalPoint recHitPos = recHit->localPosition();
0121       float res = PointExtrapolatedRPCFrame.x() - recHitPos.x();
0122       if (fabs(res) < fabs(minres)) {
0123         minres = res;
0124         cluSize = recHit->clusterSize();
0125       }
0126     }
0127     if (fabs(minres) >= (rangestrips + cluSize * 0.5) * 3) {  //3 is a typyical strip width for RPCs
0128       //std::cout<<"CSC passed, RecHits but far away "<<rpcId<<std::endl;
0129       return true;
0130     }
0131   }
0132 
0133   //std::cout<<" skiped"<<std::endl;
0134   return false;
0135 }
0136 
0137 // declare this class as a framework plugin
0138 #include "FWCore/Framework/interface/MakerMacros.h"
0139 DEFINE_FWK_MODULE(HLTRPCFilter);