CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
JetDeltaRValueMapProducer< T, C > Class Template Reference
Inheritance diagram for JetDeltaRValueMapProducer< T, C >:
edm::stream::EDProducer<>

Public Types

typedef edm::ValueMap< float > JetValueMap
 
- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Public Member Functions

 JetDeltaRValueMapProducer (edm::ParameterSet const &params)
 
 ~JetDeltaRValueMapProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

Private Member Functions

void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 

Private Attributes

const double distMax_
 
std::map< std::string, std::unique_ptr< const StringObjectFunction< C > > > evaluationMap_
 
const bool lazyParser_
 
const edm::EDGetTokenT< typename edm::View< C > > matchedToken_
 
bool multiValue_
 
const edm::EDGetTokenT< typename edm::View< T > > srcToken_
 
const std::string value_
 
const std::vector< std::string > valueLabels_
 
const std::vector< std::string > values_
 

Detailed Description

template<class T, class C = T>
class JetDeltaRValueMapProducer< T, C >

Definition at line 31 of file JetDeltaRValueMapProducer.cc.

Member Typedef Documentation

◆ JetValueMap

template<class T , class C = T>
typedef edm::ValueMap<float> JetDeltaRValueMapProducer< T, C >::JetValueMap

Definition at line 33 of file JetDeltaRValueMapProducer.cc.

Constructor & Destructor Documentation

◆ JetDeltaRValueMapProducer()

template<class T , class C = T>
JetDeltaRValueMapProducer< T, C >::JetDeltaRValueMapProducer ( edm::ParameterSet const &  params)
inline

Definition at line 35 of file JetDeltaRValueMapProducer.cc.

36  : srcToken_(consumes<typename edm::View<T> >(params.getParameter<edm::InputTag>("src"))),
37  matchedToken_(consumes<typename edm::View<C> >(params.getParameter<edm::InputTag>("matched"))),
38  distMax_(params.getParameter<double>("distMax")),
39  value_(params.existsAs<std::string>("value") ? params.getParameter<std::string>("value") : ""),
40  values_(params.existsAs<std::vector<std::string> >("values")
41  ? params.getParameter<std::vector<std::string> >("values")
42  : std::vector<std::string>()),
43  valueLabels_(params.existsAs<std::vector<std::string> >("valueLabels")
44  ? params.getParameter<std::vector<std::string> >("valueLabels")
45  : std::vector<std::string>()),
46  lazyParser_(params.existsAs<bool>("lazyParser") ? params.getParameter<bool>("lazyParser") : false),
47  multiValue_(false) {
48  if (!value_.empty()) {
49  evaluationMap_.insert(std::make_pair(
51  produces<JetValueMap>();
52  }
53 
54  if (!valueLabels_.empty() || !values_.empty()) {
55  if (valueLabels_.size() == values_.size()) {
56  multiValue_ = true;
57  for (size_t i = 0; i < valueLabels_.size(); ++i) {
58  evaluationMap_.insert(std::make_pair(
59  valueLabels_[i],
61  produces<JetValueMap>(valueLabels_[i]);
62  }
63  } else
64  edm::LogWarning("ValueLabelMismatch")
65  << "The number of value labels does not match the number of values. Values will not be evaluated.";
66  }
67  }

References JetDeltaRValueMapProducer< T, C >::evaluationMap_, mps_fire::i, JetDeltaRValueMapProducer< T, C >::lazyParser_, JetDeltaRValueMapProducer< T, C >::multiValue_, JetDeltaRValueMapProducer< T, C >::value_, JetDeltaRValueMapProducer< T, C >::valueLabels_, and JetDeltaRValueMapProducer< T, C >::values_.

◆ ~JetDeltaRValueMapProducer()

template<class T , class C = T>
JetDeltaRValueMapProducer< T, C >::~JetDeltaRValueMapProducer ( )
inlineoverride

Definition at line 69 of file JetDeltaRValueMapProducer.cc.

69 {}

Member Function Documentation

◆ produce()

template<class T , class C = T>
void JetDeltaRValueMapProducer< T, C >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
inlineoverrideprivate

Definition at line 72 of file JetDeltaRValueMapProducer.cc.

72  {
74  iEvent.getByToken(srcToken_, h_jets1);
76  iEvent.getByToken(matchedToken_, h_jets2);
77 
78  std::vector<float> values(h_jets1->size(), -99999);
79  std::map<std::string, std::vector<float> > valuesMap;
80  if (multiValue_) {
81  for (size_t i = 0; i < valueLabels_.size(); ++i)
82  valuesMap.insert(std::make_pair(valueLabels_[i], std::vector<float>(h_jets1->size(), -99999)));
83  }
84  std::vector<bool> jets1_locks(h_jets1->size(), false);
85 
86  for (typename edm::View<C>::const_iterator ibegin = h_jets2->begin(), iend = h_jets2->end(), ijet = ibegin;
87  ijet != iend;
88  ++ijet) {
89  float matched_dR2 = 1e9;
90  int matched_index = -1;
91 
92  for (typename edm::View<T>::const_iterator jbegin = h_jets1->begin(), jend = h_jets1->end(), jjet = jbegin;
93  jjet != jend;
94  ++jjet) {
95  int index = jjet - jbegin;
96 
97  if (jets1_locks.at(index))
98  continue; // skip jets that have already been matched
99 
100  float temp_dR2 = reco::deltaR2(ijet->eta(), ijet->phi(), jjet->eta(), jjet->phi());
101  if (temp_dR2 < matched_dR2) {
102  matched_dR2 = temp_dR2;
103  matched_index = index;
104  }
105  } // end loop over src jets
106 
107  if (matched_index >= 0) {
108  if (matched_dR2 > distMax_ * distMax_)
109  edm::LogInfo("MatchedJetsFarApart") << "Matched jets separated by dR greater than distMax=" << distMax_;
110  else {
111  jets1_locks.at(matched_index) = true;
112  if (!value_.empty())
113  values.at(matched_index) = (*(evaluationMap_.at(value_)))(*ijet);
114  if (multiValue_) {
115  for (size_t i = 0; i < valueLabels_.size(); ++i)
116  valuesMap.at(valueLabels_[i]).at(matched_index) = (*(evaluationMap_.at(valueLabels_[i])))(*ijet);
117  }
118  }
119  }
120  } // end loop over matched jets
121 
122  if (!value_.empty()) {
123  std::unique_ptr<JetValueMap> jetValueMap(new JetValueMap());
124 
125  JetValueMap::Filler filler(*jetValueMap);
126  filler.insert(h_jets1, values.begin(), values.end());
127  filler.fill();
128 
129  // put in Event
130  iEvent.put(std::move(jetValueMap));
131  }
132  if (multiValue_) {
133  for (size_t i = 0; i < valueLabels_.size(); ++i) {
134  std::unique_ptr<JetValueMap> jetValueMap(new JetValueMap());
135 
136  JetValueMap::Filler filler(*jetValueMap);
137  filler.insert(h_jets1, valuesMap.at(valueLabels_[i]).begin(), valuesMap.at(valueLabels_[i]).end());
138  filler.fill();
139 
140  // put in Event
141  iEvent.put(std::move(jetValueMap), valueLabels_[i]);
142  }
143  }
144  }

References reco::deltaR2(), JetDeltaRValueMapProducer< T, C >::distMax_, JetDeltaRValueMapProducer< T, C >::evaluationMap_, trigObjTnPSource_cfi::filler, mps_fire::i, iEvent, JetDeltaRValueMapProducer< T, C >::matchedToken_, eostools::move(), JetDeltaRValueMapProducer< T, C >::multiValue_, JetDeltaRValueMapProducer< T, C >::srcToken_, JetDeltaRValueMapProducer< T, C >::value_, JetDeltaRValueMapProducer< T, C >::valueLabels_, and contentValuesCheck::values.

Member Data Documentation

◆ distMax_

template<class T , class C = T>
const double JetDeltaRValueMapProducer< T, C >::distMax_
private

◆ evaluationMap_

template<class T , class C = T>
std::map<std::string, std::unique_ptr<const StringObjectFunction<C> > > JetDeltaRValueMapProducer< T, C >::evaluationMap_
private

◆ lazyParser_

template<class T , class C = T>
const bool JetDeltaRValueMapProducer< T, C >::lazyParser_
private

◆ matchedToken_

template<class T , class C = T>
const edm::EDGetTokenT<typename edm::View<C> > JetDeltaRValueMapProducer< T, C >::matchedToken_
private

◆ multiValue_

template<class T , class C = T>
bool JetDeltaRValueMapProducer< T, C >::multiValue_
private

◆ srcToken_

template<class T , class C = T>
const edm::EDGetTokenT<typename edm::View<T> > JetDeltaRValueMapProducer< T, C >::srcToken_
private

◆ value_

template<class T , class C = T>
const std::string JetDeltaRValueMapProducer< T, C >::value_
private

◆ valueLabels_

template<class T , class C = T>
const std::vector<std::string> JetDeltaRValueMapProducer< T, C >::valueLabels_
private

◆ values_

template<class T , class C = T>
const std::vector<std::string> JetDeltaRValueMapProducer< T, C >::values_
private
mps_fire.i
i
Definition: mps_fire.py:428
StringObjectFunction
Definition: StringObjectFunction.h:16
JetDeltaRValueMapProducer::distMax_
const double distMax_
Definition: JetDeltaRValueMapProducer.cc:148
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
JetDeltaRValueMapProducer::JetValueMap
edm::ValueMap< float > JetValueMap
Definition: JetDeltaRValueMapProducer.cc:33
JetDeltaRValueMapProducer::value_
const std::string value_
Definition: JetDeltaRValueMapProducer.cc:149
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
JetDeltaRValueMapProducer::valueLabels_
const std::vector< std::string > valueLabels_
Definition: JetDeltaRValueMapProducer.cc:151
edm::Handle
Definition: AssociativeIterator.h:50
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
contentValuesCheck.values
values
Definition: contentValuesCheck.py:38
JetDeltaRValueMapProducer::evaluationMap_
std::map< std::string, std::unique_ptr< const StringObjectFunction< C > > > evaluationMap_
Definition: JetDeltaRValueMapProducer.cc:154
JetDeltaRValueMapProducer::srcToken_
const edm::EDGetTokenT< typename edm::View< T > > srcToken_
Definition: JetDeltaRValueMapProducer.cc:146
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
JetDeltaRValueMapProducer::matchedToken_
const edm::EDGetTokenT< typename edm::View< C > > matchedToken_
Definition: JetDeltaRValueMapProducer.cc:147
edm::View
Definition: CaloClusterFwd.h:14
trigObjTnPSource_cfi.filler
filler
Definition: trigObjTnPSource_cfi.py:21
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
JetDeltaRValueMapProducer::values_
const std::vector< std::string > values_
Definition: JetDeltaRValueMapProducer.cc:150
iEvent
int iEvent
Definition: GenABIO.cc:224
JetDeltaRValueMapProducer::multiValue_
bool multiValue_
Definition: JetDeltaRValueMapProducer.cc:153
edm::ValueMap< float >::Filler
helper::Filler< ValueMap< float > > Filler
Definition: ValueMap.h:168
cms::cuda::device::unique_ptr
std::unique_ptr< T, impl::DeviceDeleter > unique_ptr
Definition: device_unique_ptr.h:33
JetDeltaRValueMapProducer::lazyParser_
const bool lazyParser_
Definition: JetDeltaRValueMapProducer.cc:152
eostools.move
def move(src, dest)
Definition: eostools.py:511
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::InputTag
Definition: InputTag.h:15