CMS 3D CMS Logo

OniaAddV0TracksProducer.cc
Go to the documentation of this file.
5 #include <vector>
6 
8  : events_v0{0}, total_v0{0}, total_lambda{0}, total_kshort{0} {
9  LambdaCollectionToken_ =
10  consumes<reco::VertexCompositeCandidateCollection>(ps.getParameter<edm::InputTag>("LambdaTag"));
11  KShortCollectionToken_ =
12  consumes<reco::VertexCompositeCandidateCollection>(ps.getParameter<edm::InputTag>("KShortTag"));
13 
14  produces<pat::CompositeCandidateCollection>("Kshort");
15  produces<pat::CompositeCandidateCollection>("Lambda");
16 }
17 
19  // Create unique_ptr for each collection to be stored in the Event
20  std::unique_ptr<pat::CompositeCandidateCollection> Enhanced_kShortCandidates(new pat::CompositeCandidateCollection);
21  std::unique_ptr<pat::CompositeCandidateCollection> Enhanced_lambdaCandidates(new pat::CompositeCandidateCollection);
22 
24  event.getByToken(LambdaCollectionToken_, lcandidates);
25 
27  event.getByToken(KShortCollectionToken_, kcandidates);
28 
29  int exits_l = 0;
30  int exits_k = 0;
31  for (reco::VertexCompositeCandidateCollection::const_iterator ik = kcandidates->begin(); ik != kcandidates->end();
32  ++ik) {
34  edm::RefToBase<reco::Track> ktrk0((*(dynamic_cast<const reco::RecoChargedCandidate*>(ik->daughter(0)))).track());
35  edm::RefToBase<reco::Track> ktrk1((*(dynamic_cast<const reco::RecoChargedCandidate*>(ik->daughter(1)))).track());
36  kc.addUserData<reco::Track>("track0", *ktrk0);
37  kc.addUserData<reco::Track>("track1", *ktrk1);
38  Enhanced_kShortCandidates->emplace_back(std::move(kc));
39  exits_k++;
40  }
41 
42  for (reco::VertexCompositeCandidateCollection::const_iterator il = lcandidates->begin(); il != lcandidates->end();
43  ++il) {
45  edm::RefToBase<reco::Track> ltrk0((*(dynamic_cast<const reco::RecoChargedCandidate*>(il->daughter(0)))).track());
46  edm::RefToBase<reco::Track> ltrk1((*(dynamic_cast<const reco::RecoChargedCandidate*>(il->daughter(1)))).track());
47  lc.addUserData<reco::Track>("track0", *ltrk0);
48  lc.addUserData<reco::Track>("track1", *ltrk1);
49  Enhanced_lambdaCandidates->emplace_back(std::move(lc));
50  exits_l++;
51  }
52 
53  // Write the collections to the Event
54 
55  total_v0 += exits_k;
56  total_v0 += exits_l;
57  total_kshort += exits_k;
58  total_lambda += exits_l;
59  if (exits_k || exits_l)
60  events_v0++;
61 
62  event.put(std::move(Enhanced_kShortCandidates), "Kshort");
63  event.put(std::move(Enhanced_lambdaCandidates), "Lambda");
64 }
65 
67  edm::LogVerbatim("OniaAddV0TracksSummary") << "############################\n"
68  "OniaAddV0Tracks producer report \n"
69  "############################\n"
70  "Total events with v0 : "
71  << events_v0
72  << "\n"
73  "Total v0 : "
74  << total_v0
75  << "\n"
76  "Total number of lambda : "
77  << total_lambda
78  << "\n"
79  "Total number of kshort : "
80  << total_kshort
81  << "\n"
82  "############################";
83 }
84 //define this as a plug-in
Analysis-level particle class.
Log< level::Info, true > LogVerbatim
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > KShortCollectionToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::StreamID, edm::Event &event, const edm::EventSetup &esetup) const override
std::vector< CompositeCandidate > CompositeCandidateCollection
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > LambdaCollectionToken_
OniaAddV0TracksProducer(const edm::ParameterSet &ps)
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1