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);
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_(
62  iConfig.existsAs<std::string>("objectSelection") ? iConfig.getParameter<std::string>("objectSelection") : "",
63  true) {
64  produces<edm::ValueMap<float>>();
65 }
66 
67 template <typename T>
69 
70 template <typename T>
72  using namespace edm;
73 
74  // read input
76  iEvent.getByToken(probesToken_, probes);
77 
79  iEvent.getByToken(objectsToken_, objects);
80 
81  // prepare vector for output
82  std::vector<float> values;
83 
84  // fill
85  View<reco::Candidate>::const_iterator probe, endprobes = probes->end();
86  for (probe = probes->begin(); probe != endprobes; ++probe) {
87  double dr2min = 10000;
88  for (unsigned int iObj = 0; iObj < objects->size(); iObj++) {
89  const T& obj = objects->at(iObj);
90  if (!objCut_(obj))
91  continue;
92  double dr2 = deltaR2(*probe, obj);
93  if (dr2 < dr2min) {
94  dr2min = dr2;
95  }
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 
109 // plugin definition
111 
118 
DeltaRNearestObjectComputer::objectsToken_
edm::EDGetTokenT< edm::View< T > > objectsToken_
Definition: DeltaRNearestObjectComputer.cc:72
Muon.h
configurableAnalysis::Candidate
char Candidate[]
Definition: modules.cc:20
EDProducer.h
DeltaRNearestElectronComputer
DeltaRNearestObjectComputer< reco::Electron > DeltaRNearestElectronComputer
Definition: DeltaRNearestObjectComputer.cc:113
DeltaRNearestObjectComputer::objCut_
StringCutObjectSelector< T, true > objCut_
Definition: DeltaRNearestObjectComputer.cc:73
sistrip::View
View
Definition: ConstantsForView.h:26
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
objects
Definition: __init__.py:1
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89285
Jet.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
DeltaRNearestCandidateComputer
DeltaRNearestObjectComputer< reco::Candidate > DeltaRNearestCandidateComputer
Definition: DeltaRNearestObjectComputer.cc:111
edm::Handle
Definition: AssociativeIterator.h:50
CandidateFwd.h
MakerMacros.h
Photon.h
DeltaRNearestPhotonComputer
DeltaRNearestObjectComputer< reco::Photon > DeltaRNearestPhotonComputer
Definition: DeltaRNearestObjectComputer.cc:115
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
contentValuesCheck.values
values
Definition: contentValuesCheck.py:38
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
DeltaRNearestObjectComputer::produce
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: DeltaRNearestObjectComputer.cc:70
GsfElectron.h
getGTfromDQMFile.obj
obj
Definition: getGTfromDQMFile.py:32
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
funct::true
true
Definition: Factorize.h:173
DeltaRNearestGsfElectronComputer
DeltaRNearestObjectComputer< reco::GsfElectron > DeltaRNearestGsfElectronComputer
Definition: DeltaRNearestObjectComputer.cc:114
edm::ParameterSet
Definition: ParameterSet.h:47
DeltaRNearestJetComputer
DeltaRNearestObjectComputer< reco::Jet > DeltaRNearestJetComputer
Definition: DeltaRNearestObjectComputer.cc:116
Event.h
deltaR.h
trigObjTnPSource_cfi.filler
filler
Definition: trigObjTnPSource_cfi.py:21
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup
Definition: EventSetup.h:58
hgcalPerformanceValidation.objects
objects
Definition: hgcalPerformanceValidation.py:791
DeltaRNearestObjectComputer
Definition: DeltaRNearestObjectComputer.cc:43
DeltaRNearestObjectComputer::probesToken_
edm::EDGetTokenT< edm::View< reco::Candidate > > probesToken_
Definition: DeltaRNearestObjectComputer.cc:71
DeltaRNearestMuonComputer
DeltaRNearestObjectComputer< reco::Muon > DeltaRNearestMuonComputer
Definition: DeltaRNearestObjectComputer.cc:112
InputTag.h
ValueMap.h
Electron.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
StringCutObjectSelector.h
HLTMuonOfflineAnalyzer_cfi.deltaR2
deltaR2
Definition: HLTMuonOfflineAnalyzer_cfi.py:105
T
long double T
Definition: Basic3DVectorLD.h:48
edm::ValueMap
Definition: ValueMap.h:107
StringCutObjectSelector< T, true >
edm::EDProducer
Definition: EDProducer.h:35
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
Candidate.h
View.h
ParameterSet.h
DeltaRNearestObjectComputer::DeltaRNearestObjectComputer
DeltaRNearestObjectComputer(const edm::ParameterSet &iConfig)
Definition: DeltaRNearestObjectComputer.cc:57
DeltaRNearestObjectComputer::~DeltaRNearestObjectComputer
~DeltaRNearestObjectComputer() override
Definition: DeltaRNearestObjectComputer.cc:67
edm::Event
Definition: Event.h:73