CMS 3D CMS Logo

MCTruthCompositeMatcherNew.cc
Go to the documentation of this file.
1 /* \class MCTruthCompositeMatcher
2  *
3  * \author Luca Lista, INFN
4  *
5  */
6 
17 
18 namespace reco {
19  namespace modulesNew {
20 
22  public:
23  explicit MCTruthCompositeMatcher( const edm::ParameterSet & );
24  ~MCTruthCompositeMatcher() override;
25  private:
27  std::vector<edm::EDGetTokenT<reco::GenParticleMatch> > matchMapTokens_;
28  std::vector<int> pdgId_;
29  void produce( edm::Event & , const edm::EventSetup&) override;
30  };
31 
33  srcToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("src"))),
34  matchMapTokens_( edm::vector_transform(cfg.template getParameter<std::vector<edm::InputTag> >( "matchMaps" ), [this](edm::InputTag const & tag){return consumes<reco::GenParticleMatch>(tag);} ) ),
35  pdgId_(cfg.getParameter<std::vector<int> >("matchPDGId")) {
36  produces<reco::GenParticleMatch>();
37  }
38 
40  }
41 
43  using namespace edm;
44  using namespace std;
46  evt.getByToken(srcToken_, cands);
47  size_t nMaps = matchMapTokens_.size();
48  std::vector<const GenParticleMatch *> maps;
49  maps.reserve( nMaps );
50  for( size_t i = 0; i != nMaps; ++ i ) {
52  evt.getByToken(matchMapTokens_[i], matchMap);
53  maps.push_back(& * matchMap);
54  }
56  auto matchMap = std::make_unique<GenParticleMatch>(match.ref());
57  int size = cands->size();
58  vector<int>::const_iterator begin = pdgId_.begin(), end = pdgId_.end();
59  if(size != 0) {
61  vector<int> indices(size);
62  for(int i = 0; i != size; ++ i) {
63  const Candidate & cand = (* cands)[i];
64  GenParticleRef mc = match[cand];
65  if(mc.isNull()) {
66  indices[i] = -1;
67  } else {
68  bool found = true;
69  if(begin!=end) found = find(begin, end, std::abs(mc->pdgId())) != end;
70  indices[i] = found ? int(mc.key()) : -1;
71  }
72  }
73  CandidateBaseRefProd ref(edm::makeRefToBaseProdFrom(cands->refAt(0), evt));
74  filler.insert(ref, indices.begin(), indices.end());
75  filler.fill();
76  }
77  evt.put(std::move(matchMap));
78  }
79 
80  }
81 }
82 
85 
86 namespace reco {
87  namespace modulesNew {
88 
90 
91 DEFINE_FWK_MODULE( MCTruthCompositeMatcherNew );
92 
93  }
94 }
95 
size
Write out results.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
std::vector< edm::EDGetTokenT< reco::GenParticleMatch > > matchMapTokens_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
map_type::refprod_type ref() const
reference to matched collection
key_type key() const
Accessor for product key.
Definition: Ref.h:265
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:520
#define end
Definition: vmac.h:39
bool isNull() const
Checks for null.
Definition: Ref.h:250
void produce(edm::Event &, const edm::EventSetup &) override
MCTruthCompositeMatcher MCTruthCompositeMatcherNew
RefToBaseProd< T > makeRefToBaseProdFrom(RefToBase< T > const &iRef, Event const &iEvent)
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
def move(src, dest)
Definition: eostools.py:510