CMS 3D CMS Logo

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