Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:35

0001 #include "RecoJets/JetProducers/interface/JetMuonHitsIDHelper.h"
0002 
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 
0005 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
0006 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
0007 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
0008 #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
0009 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0010 #include "DataFormats/JetReco/interface/CaloJet.h"
0011 #include "DataFormats/TrackReco/interface/Track.h"
0012 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
0013 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0014 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
0015 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
0016 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
0017 #include "DataFormats/MuonReco/interface/Muon.h"
0018 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0019 #include "DataFormats/Math/interface/deltaR.h"
0020 
0021 #include "TMath.h"
0022 #include <vector>
0023 #include <iostream>
0024 
0025 using namespace std;
0026 
0027 reco::helper::JetMuonHitsIDHelper::JetMuonHitsIDHelper(edm::ParameterSet const& pset, edm::ConsumesCollector&& iC)
0028     : trackingGeometryToken_(iC.esConsumes()) {
0029   isRECO_ = true;  // This will be "true" initially, then if the product isn't found, set to false once
0030   numberOfHits1RPC_ = 0;
0031   numberOfHits2RPC_ = 0;
0032   numberOfHits3RPC_ = 0;
0033   numberOfHits4RPC_ = 0;
0034   numberOfHitsRPC_ = 0;
0035   rpcRecHits_ = pset.getParameter<edm::InputTag>("rpcRecHits");
0036 
0037   input_rpchits_token_ = iC.consumes<RPCRecHitCollection>(rpcRecHits_);
0038 }
0039 
0040 void reco::helper::JetMuonHitsIDHelper::calculate(const edm::Event& event,
0041                                                   const edm::EventSetup& iSetup,
0042                                                   const reco::Jet& jet,
0043                                                   const int iDbg) {
0044   // initialize
0045   numberOfHits1RPC_ = 0;
0046   numberOfHits2RPC_ = 0;
0047   numberOfHits3RPC_ = 0;
0048   numberOfHits4RPC_ = 0;
0049   numberOfHitsRPC_ = 0;
0050 
0051   if (isRECO_) {  // This will be "true" initially, then if the product isn't found, set to false once
0052 
0053     // Get tracking geometry
0054     auto const& trackingGeometry = iSetup.getData(trackingGeometryToken_);
0055 
0056     //####READ RPC RecHits Collection########
0057     //#In config: RpcRecHits     = cms.InputTag("rpcRecHits")
0058     edm::Handle<RPCRecHitCollection> rpcRecHits_handle;
0059     event.getByToken(input_rpchits_token_, rpcRecHits_handle);
0060 
0061     if (!rpcRecHits_handle.isValid()) {
0062       // don't throw exception if not running on RECO
0063       edm::LogWarning("DataNotAvailable") << "JetMuonHitsIDHelper will not be run at all, this is not a RECO file.";
0064       isRECO_ = false;
0065       return;
0066     }
0067 
0068     //####calculate rpc  variables for each jet########
0069 
0070     for (RPCRecHitCollection::const_iterator itRPC = rpcRecHits_handle->begin(), itRPCEnd = rpcRecHits_handle->end();
0071          itRPC != itRPCEnd;
0072          ++itRPC) {
0073       RPCRecHit const& hit = *itRPC;
0074       DetId detid = hit.geographicalId();
0075       LocalPoint lp = hit.localPosition();
0076       const GeomDet* gd = trackingGeometry.idToDet(detid);
0077       GlobalPoint gp = gd->toGlobal(lp);
0078       double dR2 = reco::deltaR(jet.eta(), jet.phi(), static_cast<double>(gp.eta()), static_cast<double>(gp.phi()));
0079       if (dR2 < 0.5) {
0080         RPCDetId rpcChamberId = (RPCDetId)detid;
0081         numberOfHitsRPC_++;
0082         if (rpcChamberId.station() == 1)
0083           numberOfHits1RPC_++;
0084         if (rpcChamberId.station() == 2)
0085           numberOfHits2RPC_++;
0086         if (rpcChamberId.station() == 3)
0087           numberOfHits3RPC_++;
0088         if (rpcChamberId.station() == 4)
0089           numberOfHits4RPC_++;
0090       }
0091     }
0092   }
0093 }