CMS 3D CMS Logo

OniaAddV0TracksProducer.cc
Go to the documentation of this file.
5 #include <vector>
6 
8 
9  LambdaCollectionToken_ = consumes<reco::VertexCompositeCandidateCollection>(ps.getParameter<edm::InputTag>("LambdaTag"));
10  KShortCollectionToken_ = consumes<reco::VertexCompositeCandidateCollection>(ps.getParameter<edm::InputTag>("KShortTag"));
11 
12  produces< pat::CompositeCandidateCollection >("Kshort");
13  produces< pat::CompositeCandidateCollection >("Lambda");
14 
15  events_v0 = 0;
16  total_v0 = 0;
17  total_lambda = 0;
18  total_kshort = 0;
19  std::cout << "running OniaAddV0TracksProducer..." << std::endl;
20 }
21 
23 
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(); ++ik) {
38  edm::RefToBase<reco::Track> ktrk0((*(dynamic_cast<const reco::RecoChargedCandidate*>(ik->daughter(0)))).track());
39  edm::RefToBase<reco::Track> ktrk1((*(dynamic_cast<const reco::RecoChargedCandidate*>(ik->daughter(1)))).track());
40  kc->addUserData<reco::Track>( "track0", *ktrk0 );
41  kc->addUserData<reco::Track>( "track1", *ktrk1 );
42  Enhanced_kShortCandidates->push_back(*kc);
43  exits_k++;
44  }
45 
46  for (reco::VertexCompositeCandidateCollection::const_iterator il = lcandidates->begin(); il != lcandidates->end(); ++il) {
48  edm::RefToBase<reco::Track> ltrk0((*(dynamic_cast<const reco::RecoChargedCandidate*>(il->daughter(0)))).track());
49  edm::RefToBase<reco::Track> ltrk1((*(dynamic_cast<const reco::RecoChargedCandidate*>(il->daughter(1)))).track());
50  lc->addUserData<reco::Track>( "track0", *ltrk0 );
51  lc->addUserData<reco::Track>( "track1", *ltrk1 );
52  Enhanced_lambdaCandidates->push_back(*lc);
53  exits_l++;
54  }
55 
56  // Write the collections to the Event
57 
58  total_v0 += exits_k;
59  total_v0 += exits_l;
60  total_kshort += exits_k;
61  total_lambda += exits_l;
62  if (exits_k || exits_l) events_v0++;
63 
64  event.put(std::move(Enhanced_kShortCandidates),"Kshort");
65  event.put(std::move(Enhanced_lambdaCandidates),"Lambda");
66 }
67 
69  std::cout << "############################" << std::endl;
70  std::cout << "OniaAddV0Tracks producer report " << std::endl;
71  std::cout << "############################" << std::endl;
72  std::cout << "Total events with v0 : " << events_v0 << std::endl;
73  std::cout << "Total v0 : " << total_v0 << std::endl;
74  std::cout << "Total number of lambda : " << total_lambda << std::endl;
75  std::cout << "Total number of kshort : " << total_kshort << std::endl;
76  std::cout << "############################" << std::endl;
77 }
78 //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:17
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:309
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1