CMS 3D CMS Logo

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  std::vector<edm::InputTag> inputCollVector_ = iConfig.getParameter<std::vector<InputTag>>("InputCollections");
9  for (unsigned int ii = 0; ii < inputCollVector_.size(); ++ii) {
10  inputColl_.push_back(consumes<LorentzVectorCollection>(inputCollVector_[ii]));
11  }
12  matchDeltaR_ = iConfig.getParameter<double>("MatchDeltaR");
13  outName_ = iConfig.getParameter<string>("OutputCollection");
14 
15  produces<LorentzVectorCollection>(outName_);
16 }
17 
19 
21  unique_ptr<LorentzVectorCollection> out_product(new LorentzVectorCollection);
22 
23  // Create The Handles..
24  std::vector<Handle<LorentzVectorCollection>> handles;
25 
26  bool allCollectionsExist = true;
27  // Map the Handles to the collections if all collections exist
28  for (size_t i = 0; i < inputColl_.size(); ++i) {
30  if (iEvent.getByToken(inputColl_[i], tmp)) {
31  handles.push_back(tmp);
32  } else {
33  allCollectionsExist = false;
34  }
35  }
36 
37  // The reference output object collection will be the first one..
38  if (allCollectionsExist) {
39  // loop on the first collection
40  for (size_t i = 0; i < (handles[0])->size(); ++i) {
41  bool MatchedObj = true;
42 
43  // get reference Vector
44  const LorentzVector lvRef = (*(handles[0]))[i];
45 
46  // Loop on all other collections and match
47  for (size_t j = 1; j < handles.size(); ++j) {
48  if (!match(lvRef, *(handles[j])))
49  MatchedObj = false;
50  }
51 
52  // If the Object is Matched Everywhere store
53  if (MatchedObj) {
54  out_product->push_back(lvRef);
55  }
56  }
57 
58  // Put product to file
59  iEvent.put(std::move(out_product), outName_);
60  }
61 }
62 
64  bool matched = false;
65 
66  if (!lvcol.empty())
67  for (LorentzVectorCollection::const_iterator it = lvcol.begin(); it != lvcol.end(); ++it) {
68  double delta = ROOT::Math::VectorUtil::DeltaR(lv, *it);
69  if (delta < matchDeltaR_) {
70  matched = true;
71  }
72  }
73 
74  return matched;
75 }
size
Write out results.
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool match(const LorentzVector &, const LorentzVectorCollection &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< LorentzVector > LorentzVectorCollection
int iEvent
Definition: GenABIO.cc:224
ii
Definition: cuy.py:590
HLTTauRefCombiner(const edm::ParameterSet &)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
HLT enums.
~HLTTauRefCombiner() override
math::XYZTLorentzVectorD LorentzVector
def move(src, dest)
Definition: eostools.py:511
std::string match(BranchDescription const &a, BranchDescription const &b, std::string const &fileName)
void produce(edm::Event &, const edm::EventSetup &) override