CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/special/src/HLTRPCFilter.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Class:      HLTRPCFilter
00004 // 
00012 //
00013 // Original Author:  Camilo Andres Carrillo Montoya
00014 //         Created:  Thu Oct 29 11:04:22 CET 2009
00015 // $Id: HLTRPCFilter.cc,v 1.2 2010/02/23 10:20:56 carrillo Exp $
00016 //
00017 //
00018 
00019 #include "HLTrigger/special/interface/HLTRPCFilter.h"
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 //
00025 // constants, enums and typedefs
00026 //
00027 
00028 //
00029 // static data member definitions
00030 //
00031 
00032 //
00033 // constructors and destructor
00034 //
00035 
00036 HLTRPCFilter::HLTRPCFilter(const edm::ParameterSet& iConfig)
00037 {
00038    //now do what ever initialization is needed
00039 
00040   rangestrips = iConfig.getUntrackedParameter<double>("rangestrips",1.);
00041   rpcRecHitsLabel = iConfig.getParameter<edm::InputTag>("rpcRecHits");
00042   rpcDTPointsLabel  = iConfig.getParameter<edm::InputTag>("rpcDTPoints");
00043   rpcCSCPointsLabel  = iConfig.getParameter<edm::InputTag>("rpcCSCPoints");
00044   
00045 }
00046 
00047 
00048 HLTRPCFilter::~HLTRPCFilter()
00049 {
00050  
00051    // do anything here that needs to be done at desctruction time
00052    // (e.g. close files, deallocate resources etc.)
00053 
00054 }
00055 
00056 
00057 //
00058 // member functions
00059 //
00060 
00061 // ------------ method called on each new Event  ------------
00062 bool HLTRPCFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00063 {
00064   edm::Handle<RPCRecHitCollection> rpcHits;
00065   iEvent.getByLabel(rpcRecHitsLabel,rpcHits);
00066 
00067   RPCRecHitCollection::const_iterator rpcPoint;
00068  
00069   if(rpcHits->begin()==rpcHits->end()){
00070     //std::cout<<" skiped preventing no RPC runs"<<std::endl;
00071     return false;
00072   }
00073 
00074   edm::Handle<RPCRecHitCollection> rpcDTPoints;
00075   iEvent.getByLabel(rpcDTPointsLabel,rpcDTPoints);
00076 
00077   edm::Handle<RPCRecHitCollection> rpcCSCPoints;
00078   iEvent.getByLabel(rpcCSCPointsLabel,rpcCSCPoints);
00079 
00080   float cluSize = 0;
00081   
00082   //DTPart
00083 
00084   for(rpcPoint = rpcDTPoints->begin(); rpcPoint != rpcDTPoints->end(); rpcPoint++){
00085     LocalPoint PointExtrapolatedRPCFrame = rpcPoint->localPosition();
00086     RPCDetId  rpcId = rpcPoint->rpcId();
00087     typedef std::pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
00088     rangeRecHits recHitCollection =  rpcHits->get(rpcId);
00089     if(recHitCollection.first==recHitCollection.second){
00090       //std::cout<<"DT passed, no rechits for this RPCId =  "<<rpcId<<std::endl;
00091       return true;
00092     }
00093     float minres = 3000.;
00094     RPCRecHitCollection::const_iterator recHit;
00095     for(recHit = recHitCollection.first; recHit != recHitCollection.second ; recHit++) {
00096       LocalPoint recHitPos=recHit->localPosition();
00097       float res=PointExtrapolatedRPCFrame.x()- recHitPos.x();
00098       if(fabs(res)<fabs(minres)){
00099         minres=res;
00100         cluSize = recHit->clusterSize();
00101       }
00102     }
00103     if(fabs(minres)>=(rangestrips+cluSize*0.5)*3){ //3 is a typyical strip width for RPCs
00104       //std::cout<<"DT passed, RecHits but far away "<<rpcId<<std::endl;
00105       return true;
00106     }
00107   }
00108 
00109   //CSCPart
00110 
00111   for(rpcPoint = rpcCSCPoints->begin(); rpcPoint != rpcCSCPoints->end(); rpcPoint++){
00112     LocalPoint PointExtrapolatedRPCFrame = rpcPoint->localPosition();
00113     RPCDetId  rpcId = rpcPoint->rpcId();
00114     typedef std::pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
00115     rangeRecHits recHitCollection =  rpcHits->get(rpcId);
00116     if(recHitCollection.first==recHitCollection.second){
00117       //std::cout<<"CSC passed, no rechits for this RPCId =  "<<rpcId<<std::endl;
00118       return true;
00119     }
00120     float minres = 3000.;
00121     RPCRecHitCollection::const_iterator recHit;
00122     for(recHit = recHitCollection.first; recHit != recHitCollection.second ; recHit++) {
00123       LocalPoint recHitPos=recHit->localPosition();
00124       float res=PointExtrapolatedRPCFrame.x()- recHitPos.x();
00125       if(fabs(res)<fabs(minres)){
00126         minres=res;
00127         cluSize = recHit->clusterSize();
00128       }
00129     }
00130     if(fabs(minres)>=(rangestrips+cluSize*0.5)*3){ //3 is a typyical strip width for RPCs
00131       //std::cout<<"CSC passed, RecHits but far away "<<rpcId<<std::endl;
00132       return true;
00133     }
00134   }
00135 
00136   //std::cout<<" skiped"<<std::endl;
00137   return false;
00138 }
00139 
00140 
00141 
00142 
00143 // ------------ method called once each job just before starting event loop  ------------
00144 void 
00145 HLTRPCFilter::beginJob()
00146 {
00147 }
00148 
00149 // ------------ method called once each job just after ending the event loop  ------------
00150 void 
00151 HLTRPCFilter::endJob() {
00152 }
00153 
00154 //define this as a plug-in
00155 DEFINE_FWK_MODULE(HLTRPCFilter);