CMS 3D CMS Logo

DeltaRNearestObjectComputer.cc
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: CMS detector at the CERN
3  *
4  * Package: PhysicsTools/TagAndProbe
5  *
6  *
7  * Authors:
8  * Giovanni Petrucciani, UCSD - Giovanni.Petrucciani@cern.ch
9  *
10  * Description:
11  * - Matches a given object with other objects using deltaR-matching.
12  * - For example: can match a photon with track within a given deltaR.
13  * - Saves collection of the reference vectors of matched objects.
14  * History:
15  *
16  * Kalanand Mishra, Fermilab - kalanand@fnal.gov
17  * Extended the class to compute deltaR with respect to any object
18  * (i.e., Candidate, Jet, Muon, Electron, or Photon). The previous
19  * version of this class could deltaR only with respect to reco::Candidates.
20  * This didn't work if one wanted to apply selection cuts on the Candidate's
21  * RefToBase object.
22  *****************************************************************************/
23 
28 
32 
40 
42 
43 template<typename T>
45  public:
46  explicit DeltaRNearestObjectComputer(const edm::ParameterSet & iConfig);
47  ~DeltaRNearestObjectComputer() override ;
48 
49  void produce(edm::Event & iEvent, const edm::EventSetup& iSetup) override;
50 
51  private:
54  StringCutObjectSelector<T,true> objCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
55 };
56 
57 template<typename T>
59  probesToken_(consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("probes"))),
60  objectsToken_(consumes<edm::View<T> >(iConfig.getParameter<edm::InputTag>("objects"))),
61  objCut_(iConfig.existsAs<std::string>("objectSelection") ? iConfig.getParameter<std::string>("objectSelection") : "", true)
62 {
63  produces<edm::ValueMap<float> >();
64 }
65 
66 
67 template<typename T>
69 {
70 }
71 
72 template<typename T>
73 void
75  using namespace edm;
76 
77  // read input
79  iEvent.getByToken(probesToken_, probes);
80 
82  iEvent.getByToken(objectsToken_, objects);
83 
84  // prepare vector for output
85  std::vector<float> values;
86 
87  // fill
88  View<reco::Candidate>::const_iterator probe, endprobes = probes->end();
89  for (probe = probes->begin(); probe != endprobes; ++probe) {
90  double dr2min = 10000;
91  for (unsigned int iObj=0; iObj<objects->size(); iObj++) {
92  const T& obj = objects->at(iObj);
93  if (!objCut_(obj)) continue;
94  double dr2 = deltaR2(*probe, obj);
95  if (dr2 < dr2min) { dr2min = dr2; }
96  }
97  values.push_back(sqrt(dr2min));
98  }
99 
100  // convert into ValueMap and store
101  auto valMap = std::make_unique<ValueMap<float>>();
103  filler.insert(probes, values.begin(), values.end());
104  filler.fill();
105  iEvent.put(std::move(valMap));
106 }
107 
108 
110 // plugin definition
112 
119 
DeltaRNearestObjectComputer< reco::GsfElectron > DeltaRNearestGsfElectronComputer
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
StringCutObjectSelector< T, true > objCut_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
DeltaRNearestObjectComputer< reco::Muon > DeltaRNearestMuonComputer
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edm::View< T > > objectsToken_
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:18
DeltaRNearestObjectComputer< reco::Candidate > DeltaRNearestCandidateComputer
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
DeltaRNearestObjectComputer< reco::Electron > DeltaRNearestElectronComputer
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
DeltaRNearestObjectComputer< reco::Jet > DeltaRNearestJetComputer
DeltaRNearestObjectComputer(const edm::ParameterSet &iConfig)
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
edm::EDGetTokenT< edm::View< reco::Candidate > > probesToken_
long double T
DeltaRNearestObjectComputer< reco::Photon > DeltaRNearestPhotonComputer
def move(src, dest)
Definition: eostools.py:511