Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Class:      HLTRPCTrigNoSyncFilter
0004 //
0005 /**\class RPCPathChambFilter HLTRPCTrigNoSyncFilter.cc HLTRPCTrigNoSyncFilter.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 "HLTRPCTrigNoSyncFilter.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 typedef struct {
0036   int id;
0037   int bx;
0038   GlobalPoint gp;
0039 } RPC4DHit;
0040 
0041 bool bigmag(const RPC4DHit& Point1, const RPC4DHit& Point2) {
0042   if ((Point2).gp.mag() > (Point1).gp.mag())
0043     return true;
0044   else
0045     return false;
0046 }
0047 
0048 HLTRPCTrigNoSyncFilter::HLTRPCTrigNoSyncFilter(const edm::ParameterSet& iConfig)
0049     : HLTFilter(iConfig), muonGeometryRecordToken_(esConsumes()) {
0050   //now do what ever initialization is needed
0051   m_GMTInputTag = iConfig.getParameter<edm::InputTag>("GMTInputTag");
0052   rpcRecHitsLabel = iConfig.getParameter<edm::InputTag>("rpcRecHits");
0053   m_GMTInputToken = consumes<L1MuGMTReadoutCollection>(m_GMTInputTag);
0054   rpcRecHitsToken = consumes<RPCRecHitCollection>(rpcRecHitsLabel);
0055 }
0056 
0057 HLTRPCTrigNoSyncFilter::~HLTRPCTrigNoSyncFilter() {
0058   // do anything here that needs to be done at desctruction time
0059   // (e.g. close files, deallocate resources etc.)
0060 }
0061 
0062 void HLTRPCTrigNoSyncFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0063   edm::ParameterSetDescription desc;
0064   makeHLTFilterDescription(desc);
0065   desc.add<edm::InputTag>("GMTInputTag", edm::InputTag("hltGtDigis"));
0066   desc.add<edm::InputTag>("rpcRecHits", edm::InputTag("hltRpcRecHits"));
0067   descriptions.add("hltRPCTrigNoSyncFilter", desc);
0068 }
0069 
0070 //
0071 // member functions
0072 //
0073 
0074 // ------------ method called on each new Event  ------------
0075 bool HLTRPCTrigNoSyncFilter::hltFilter(edm::Event& iEvent,
0076                                        const edm::EventSetup& iSetup,
0077                                        trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0078   std::vector<RPC4DHit> GlobalRPC4DHits;
0079   std::vector<RPC4DHit> GlobalRPC4DHitsNoBx0;
0080 
0081   edm::Handle<RPCRecHitCollection> rpcRecHits;
0082 
0083   //std::cout <<"\t Getting the RPC RecHits"<<std::endl;
0084 
0085   iEvent.getByToken(rpcRecHitsToken, rpcRecHits);
0086 
0087   if (!rpcRecHits.isValid()) {
0088     //std::cout<<"no valid RPC Collection"<<std::endl;
0089     //std::cout<<"event skipped"<<std::endl;
0090     return false;
0091   }
0092   if (rpcRecHits->begin() == rpcRecHits->end()) {
0093     //std::cout<<"no RPC Hits in this event"<<std::endl;
0094     //std::cout<<"event skipped"<<std::endl;
0095     return false;
0096   }
0097 
0098   RPCRecHitCollection::const_iterator recHit;
0099 
0100   auto const& rpcGeo = iSetup.getHandle(muonGeometryRecordToken_);
0101 
0102   int k = 0;
0103 
0104   for (recHit = rpcRecHits->begin(); recHit != rpcRecHits->end(); recHit++) {
0105     RPCDetId rollId = (RPCDetId)(*recHit).rpcId();
0106     LocalPoint recHitPos = recHit->localPosition();
0107     const RPCRoll* rollasociated = rpcGeo->roll(rollId);
0108     const BoundPlane& RPCSurface = rollasociated->surface();
0109     GlobalPoint RecHitInGlobal = RPCSurface.toGlobal(recHitPos);
0110 
0111     int BX = recHit->BunchX();
0112     //std::cout<<"\t \t We have an RPC Rec Hit! bx="<<BX<<" Roll="<<rollId<<" Global Position="<<RecHitInGlobal<<std::endl;
0113 
0114     RPC4DHit ThisHit;
0115     ThisHit.bx = BX;
0116     ThisHit.gp = RecHitInGlobal;
0117     ThisHit.id = k;
0118     GlobalRPC4DHits.push_back(ThisHit);
0119     if (BX != 0)
0120       GlobalRPC4DHitsNoBx0.push_back(ThisHit);
0121     k++;
0122   }
0123 
0124   if (GlobalRPC4DHitsNoBx0.empty()) {
0125     //std::cout<<"all RPC Hits are syncrhonized"<<std::endl;
0126     //std::cout<<"event skipped"<<std::endl;
0127     return false;
0128   }
0129 
0130   if (GlobalRPC4DHitsNoBx0.size() > 100) {
0131     //std::cout<<"too many rpcHits preventing HLT eternal loop"<<std::endl;
0132     //std::cout<<"event skipped"<<std::endl;
0133     return false;
0134   }
0135 
0136   edm::Handle<L1MuGMTReadoutCollection> gmtrc_handle;
0137   iEvent.getByToken(m_GMTInputToken, gmtrc_handle);
0138 
0139   std::vector<L1MuGMTExtendedCand> gmt_candidates = (*gmtrc_handle).getRecord().getGMTCands();
0140 
0141   std::vector<L1MuGMTExtendedCand>::const_iterator candidate;
0142   //std::cout<<"The number of GMT Candidates in this event is "<<gmt_candidates.size()<<std::endl;
0143 
0144   if (gmt_candidates.empty()) {
0145     //std::cout<<"event skipped no gmt candidates"<<std::endl;
0146     return false;
0147   }
0148 
0149   for (candidate = gmt_candidates.begin(); candidate != gmt_candidates.end(); candidate++) {
0150     int qual = candidate->quality();
0151     //std::cout<<"quality of this GMT Candidate (should be >5)= "<<qual<<std::endl;
0152     if (qual < 5)
0153       continue;
0154     float eta = candidate->etaValue();
0155     float phi = candidate->phiValue();
0156 
0157     //std::cout<<" phi="<<phi<<" eta="<<eta<<std::endl;
0158 
0159     //Creating container in this etaphi direction;
0160 
0161     std::vector<RPC4DHit> PointsForGMT;
0162 
0163     for (auto& Point : GlobalRPC4DHitsNoBx0) {
0164       float phiP = Point.gp.phi();
0165       float etaP = Point.gp.eta();
0166       //std::cout<<"\t \t GMT   phi="<<phi<<" eta="<<eta<<std::endl;
0167       //std::cout<<"\t \t Point phi="<<phiP<<" eta="<< etaP<<std::endl;
0168       //std::cout<<"\t \t "<<fabs(phi-phiP)<<" < 0,1? "<<fabs(eta-etaP)<<" < 0.20 ?"<<std::endl;
0169       if (fabs(phi - phiP) <= 0.1 && fabs(eta - etaP) <= 0.20) {
0170         PointsForGMT.push_back(Point);
0171         //std::cout<<"\t \t match!"<<std::endl;
0172       }
0173     }
0174 
0175     //std::cout<<"\t Points matching this GMT="<<PointsForGMT.size()<<std::endl;
0176 
0177     if (PointsForGMT.empty()) {
0178       //std::cout<<"\t Not enough RPCRecHits not syncrhonized for this GMT!!!"<<std::endl;
0179       continue;
0180     }
0181 
0182     std::sort(PointsForGMT.begin(), PointsForGMT.end(), bigmag);
0183 
0184     //std::cout<<"GMT candidate phi="<<phi<<std::endl;
0185     //std::cout<<"GMT candidate eta="<<eta<<std::endl;
0186 
0187     int lastbx = -7;
0188     bool outOfTime = false;
0189     bool incr = true;
0190 
0191     //std::cout<<"\t \t loop on the RPCHit4D!!!"<<std::endl;
0192     for (auto point = PointsForGMT.begin(); point < PointsForGMT.end(); ++point) {
0193       //float r=point->gp.mag();
0194       outOfTime |= (point->bx != 0);  //condition 1: at least one measurement must have BX!=0
0195       incr &= (point->bx >= lastbx);  //condition 2: BX must be increase when going inside-out.
0196       lastbx = point->bx;
0197     }
0198     //std::cout<<"\t \t";
0199 
0200     //for(std::vector<RPC4DHit>::iterator point = PointsForGMT.begin(); point < PointsForGMT.end(); ++point) {
0201     //std::cout<<point->bx;
0202     //}
0203     //std::cout<<std::endl;
0204 
0205     bool Candidate = (outOfTime && incr);
0206 
0207     if (Candidate) {
0208       //std::cout<<" Event Passed We found an strange trigger Candidate phi="<<phi<<" eta="<<eta<<std::endl;
0209       return true;
0210     }
0211   }
0212 
0213   //std::cout<<"event skipped rechits out of time but not related with a GMT "<<std::endl;
0214   return false;
0215 }
0216 
0217 // ------------ method called once each job just before starting event loop  ------------
0218 void HLTRPCTrigNoSyncFilter::beginJob() {}
0219 
0220 // ------------ method called once each job just after ending the event loop  ------------
0221 void HLTRPCTrigNoSyncFilter::endJob() {}
0222 
0223 // declare this class as a framework plugin
0224 DEFINE_FWK_MODULE(HLTRPCTrigNoSyncFilter);