CMS 3D CMS Logo

HLTTauRefCombiner.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
3 
4 
5 
6 using namespace edm;
7 using namespace std;
8 
10 {
11  std::vector<edm::InputTag> inputCollVector_ = iConfig.getParameter< std::vector<InputTag> >("InputCollections");
12  for(unsigned int ii=0; ii<inputCollVector_.size(); ++ii)
13  {
14  inputColl_.push_back( consumes<LorentzVectorCollection>(inputCollVector_[ii]) );
15  }
16  matchDeltaR_ = iConfig.getParameter<double>("MatchDeltaR");
17  outName_ = iConfig.getParameter<string>("OutputCollection");
18 
19  produces<LorentzVectorCollection> ( outName_);
20 }
21 
23 
25 {
26  unique_ptr<LorentzVectorCollection> out_product(new LorentzVectorCollection);
27 
28  //Create The Handles..
29  std::vector< Handle<LorentzVectorCollection> > handles;
30 
31  bool allCollectionsExist = true;
32  //Map the Handles to the collections if all collections exist
33  for(size_t i = 0;i<inputColl_.size();++i)
34  {
36  if(iEvent.getByToken(inputColl_[i],tmp))
37  {
38  handles.push_back(tmp);
39  }
40  else
41  {
42  allCollectionsExist = false;
43  }
44 
45  }
46 
47  //The reference output object collection will be the first one..
48  if(allCollectionsExist)
49  {
50  //loop on the first collection
51  for(size_t i = 0; i < (handles[0])->size();++i)
52  {
53  bool MatchedObj = true;
54 
55  //get reference Vector
56  const LorentzVector lvRef = (*(handles[0]))[i];
57 
58  //Loop on all other collections and match
59  for(size_t j = 1; j < handles.size();++j)
60  {
61  if(!match(lvRef,*(handles[j])))
62  MatchedObj = false;
63 
64  }
65 
66  //If the Object is Matched Everywhere store
67  if(MatchedObj)
68  {
69  out_product->push_back(lvRef);
70  }
71 
72 
73  }
74 
75  //Put product to file
76  iEvent.put(std::move(out_product),outName_);
77 
78  }
79 
80 
81 
82 
83 }
84 
85 
86 
87 bool
89 {
90  bool matched=false;
91 
92  if(!lvcol.empty())
93  for(LorentzVectorCollection::const_iterator it = lvcol.begin();it!=lvcol.end();++it)
94  {
95  double delta = ROOT::Math::VectorUtil::DeltaR(lv,*it);
96  if(delta<matchDeltaR_)
97  {
98  matched=true;
99 
100  }
101  }
102 
103 
104 
105  return matched;
106 }
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:137
bool match(const LorentzVector &, const LorentzVectorCollection &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
std::vector< LorentzVector > LorentzVectorCollection
int iEvent
Definition: GenABIO.cc:230
ii
Definition: cuy.py:589
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:510
std::string match(BranchDescription const &a, BranchDescription const &b, std::string const &fileName)
void produce(edm::Event &, const edm::EventSetup &) override