CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTTauRefCombiner.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
3 
4 using namespace edm;
5 using namespace std;
6 
8  : matchDeltaR_{iConfig.getParameter<double>("MatchDeltaR")},
9  outName_{iConfig.getParameter<string>("OutputCollection")} {
10  std::vector<edm::InputTag> inputCollVector_ = iConfig.getParameter<std::vector<InputTag>>("InputCollections");
11  for (unsigned int ii = 0; ii < inputCollVector_.size(); ++ii) {
12  inputColl_.push_back(consumes<LorentzVectorCollection>(inputCollVector_[ii]));
13  }
14 
15  produces<LorentzVectorCollection>(outName_);
16 }
17 
19  unique_ptr<LorentzVectorCollection> out_product(new LorentzVectorCollection);
20 
21  // Create The Handles..
22  std::vector<Handle<LorentzVectorCollection>> handles;
23 
24  bool allCollectionsExist = true;
25  // Map the Handles to the collections if all collections exist
26  for (size_t i = 0; i < inputColl_.size(); ++i) {
28  if (iEvent.getByToken(inputColl_[i], tmp)) {
29  handles.push_back(tmp);
30  } else {
31  allCollectionsExist = false;
32  }
33  }
34 
35  // The reference output object collection will be the first one..
36  if (allCollectionsExist) {
37  // loop on the first collection
38  for (size_t i = 0; i < (handles[0])->size(); ++i) {
39  bool MatchedObj = true;
40 
41  // get reference Vector
42  const LorentzVector lvRef = (*(handles[0]))[i];
43 
44  // Loop on all other collections and match
45  for (size_t j = 1; j < handles.size(); ++j) {
46  if (!match(lvRef, *(handles[j])))
47  MatchedObj = false;
48  }
49 
50  // If the Object is Matched Everywhere store
51  if (MatchedObj) {
52  out_product->push_back(lvRef);
53  }
54  }
55 
56  // Put product to file
57  iEvent.put(std::move(out_product), outName_);
58  }
59 }
60 
61 bool HLTTauRefCombiner::match(const LorentzVector &lv, const LorentzVectorCollection &lvcol) const {
62  bool matched = false;
63 
64  if (!lvcol.empty())
65  for (LorentzVectorCollection::const_iterator it = lvcol.begin(); it != lvcol.end(); ++it) {
66  double delta = ROOT::Math::VectorUtil::DeltaR(lv, *it);
67  if (delta < matchDeltaR_) {
68  matched = true;
69  }
70  }
71 
72  return matched;
73 }
const double matchDeltaR_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const std::string outName_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
std::vector< LorentzVector > LorentzVectorCollection
int ii
Definition: cuy.py:589
bool match(const LorentzVector &, const LorentzVectorCollection &) const
int iEvent
Definition: GenABIO.cc:224
std::vector< BaseVolumeHandle * > handles
def move
Definition: eostools.py:511
HLTTauRefCombiner(const edm::ParameterSet &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< edm::EDGetTokenT< LorentzVectorCollection > > inputColl_
math::XYZTLorentzVectorD LorentzVector
tmp
align.sh
Definition: createJobs.py:716
tuple size
Write out results.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override