CMS 3D CMS Logo

GenPUProtonProducer.cc
Go to the documentation of this file.
1 /* \class GenPUProtonProducer
2  *
3  * Modification of GenParticleProducer.
4  * Saves final state protons from HepMC events in Crossing Frame, in the generator-particle format.
5  *
6  * Note: Use the option USER_CXXFLAGS=-DEDM_ML_DEBUG with SCRAM in order to enable debug messages.
7  *
8  * March 9, 2017 : Initial version.
9  * March 14, 2017 : Updated debug messages.
10  * July 27, 2017 : Removed extra loop initially inherited from GenParticleProducer.
11  * August 17, 2017 : Replaced std::auto_ptr with std::unique_ptr.
12  * September 6, 2017 : Updated module to edm::global::EDProducer with ConvertParticle as RunCache following GenParticleProducer.
13  *
14  */
15 
20 
21 #include <vector>
22 #include <string>
23 #include <unordered_map>
24 
25 namespace {
26 
27  class ConvertParticle {
28  public:
29  static constexpr int PDGCacheMax = 32768;
30  static constexpr double mmToCm = 0.1;
31 
32  ConvertParticle()
33  : abortOnUnknownPDGCode_(true), initialized_(false), chargeP_(PDGCacheMax, 0), chargeM_(PDGCacheMax, 0){};
34 
35  ConvertParticle(bool abortOnUnknownPDGCode)
36  : abortOnUnknownPDGCode_(abortOnUnknownPDGCode),
37  initialized_(false),
38  chargeP_(PDGCacheMax, 0),
39  chargeM_(PDGCacheMax, 0){};
40 
41  ~ConvertParticle(){};
42 
43  bool initialized() const { return initialized_; }
44 
45  void init(HepPDT::ParticleDataTable const& pdt) {
46  if (!initialized_) {
47  for (HepPDT::ParticleDataTable::const_iterator p = pdt.begin(); p != pdt.end(); ++p) {
48  HepPDT::ParticleID const& id = p->first;
49  int pdgId = id.pid(), apdgId = std::abs(pdgId);
50  int q3 = id.threeCharge();
51  if (apdgId < PDGCacheMax && pdgId > 0) {
52  chargeP_[apdgId] = q3;
53  chargeM_[apdgId] = -q3;
54  } else if (apdgId < PDGCacheMax) {
55  chargeP_[apdgId] = -q3;
56  chargeM_[apdgId] = q3;
57  } else {
58  chargeMap_.emplace(pdgId, q3);
59  chargeMap_.emplace(-pdgId, -q3);
60  }
61  }
62  initialized_ = true;
63  }
64  }
65 
66  bool operator()(reco::GenParticle& cand, HepMC::GenParticle const* part) const {
67  reco::Candidate::LorentzVector p4(part->momentum());
68  int pdgId = part->pdg_id();
69  cand.setThreeCharge(chargeTimesThree(pdgId));
70  cand.setPdgId(pdgId);
71  cand.setStatus(part->status());
72  cand.setP4(p4);
73  cand.setCollisionId(0);
74  HepMC::GenVertex const* v = part->production_vertex();
75  if (v != nullptr) {
76  HepMC::ThreeVector vtx = v->point3d();
77  reco::Candidate::Point vertex(vtx.x() * mmToCm, vtx.y() * mmToCm, vtx.z() * mmToCm);
78  cand.setVertex(vertex);
79  } else {
80  cand.setVertex(reco::Candidate::Point(0, 0, 0));
81  }
82  return true;
83  }
84 
85  private:
86  bool abortOnUnknownPDGCode_;
87  bool initialized_;
88  std::vector<int> chargeP_, chargeM_;
89  std::unordered_map<int, int> chargeMap_;
90 
91  int chargeTimesThree(int id) const {
92  if (std::abs(id) < PDGCacheMax)
93  return id > 0 ? chargeP_[id] : chargeM_[-id];
94 
95  auto f = chargeMap_.find(id);
96  if (f == chargeMap_.end()) {
97  if (abortOnUnknownPDGCode_)
98  throw edm::Exception(edm::errors::LogicError) << "invalid PDG id: " << id << std::endl;
99  else
100  return HepPDT::ParticleID(id).threeCharge();
101  }
102  return f->second;
103  }
104  };
105 
106  class SelectProton {
107  public:
108  bool operator()(HepMC::GenParticle const* part, double minPz) const {
109  bool selection = ((!part->end_vertex() && part->status() == 1) && (part->pdg_id() == 2212) &&
110  (TMath::Abs(part->momentum().pz()) >= minPz));
111  return selection;
112  }
113  };
114 
115 } // Anonymous namespace
116 
121 
127 
128 namespace edm {
129  class ParameterSet;
130 }
131 
132 class GenPUProtonProducer : public edm::global::EDProducer<edm::RunCache<ConvertParticle> > {
133 public:
135  ~GenPUProtonProducer() override;
136 
137  void produce(edm::StreamID, edm::Event& e, const edm::EventSetup&) const override;
138  std::shared_ptr<ConvertParticle> globalBeginRun(const edm::Run&, const edm::EventSetup&) const override;
139  void globalEndRun(edm::Run const&, edm::EventSetup const&) const override{};
140 
141 private:
143 
145  std::vector<int> bunchList_;
146  double minPz_;
147  SelectProton select_;
148 };
149 
157 
161 
162 #include <algorithm>
163 
164 using namespace edm;
165 using namespace reco;
166 using namespace std;
167 using namespace HepMC;
168 
170  : abortOnUnknownPDGCode_(cfg.getUntrackedParameter<bool>("abortOnUnknownPDGCode", true)),
171  bunchList_(cfg.getParameter<vector<int> >("bunchCrossingList")),
172  minPz_(cfg.getParameter<double>("minPz")) {
173  produces<GenParticleCollection>();
174  mixToken_ =
175  consumes<CrossingFrame<HepMCProduct> >(InputTag(cfg.getParameter<std::string>("mix"), "generatorSmeared"));
176 }
177 
179 
180 std::shared_ptr<ConvertParticle> GenPUProtonProducer::globalBeginRun(const Run&, const EventSetup& es) const {
182  es.getData(pdt);
183  auto convert_ptr = std::make_shared<ConvertParticle>(abortOnUnknownPDGCode_);
184  if (!convert_ptr->initialized())
185  convert_ptr->init(*pdt);
186 
187  return convert_ptr;
188 }
189 
191  size_t totalSize = 0;
192  size_t npiles = 1;
193 
195  evt.getByToken(mixToken_, cf);
196  std::unique_ptr<MixCollection<HepMCProduct> > cfhepmcprod(new MixCollection<HepMCProduct>(cf.product()));
197  npiles = cfhepmcprod->size();
198 
199  LogDebug("GenPUProtonProducer") << " Number of pile-up events : " << npiles << endl;
200 
201  for (size_t icf = 0; icf < npiles; ++icf) {
202  LogDebug("GenPUProtonProducer") << "CF " << icf
203  << " size : " << cfhepmcprod->getObject(icf).GetEvent()->particles_size() << endl;
204  totalSize += cfhepmcprod->getObject(icf).GetEvent()->particles_size();
205  }
206  LogDebug("GenPUProtonProducer") << "Total size : " << totalSize << endl;
207 
208  // Initialise containers
209  auto candsPtr = std::make_unique<GenParticleCollection>();
210  GenParticleCollection& cands = *candsPtr;
211 
212  // Loop over pile-up events
213  ConvertParticle const& convertParticle_ = *runCache(evt.getRun().index());
214 
216  unsigned int total_number_of_protons = 0;
217  size_t idx_mix = 0;
218  // Fill collection
219  for (mixHepMC_itr = cfhepmcprod->begin(); mixHepMC_itr != cfhepmcprod->end(); ++mixHepMC_itr, ++idx_mix) {
220  int bunch = mixHepMC_itr.bunch();
221 
222  if (find(bunchList_.begin(), bunchList_.end(), bunch) != bunchList_.end()) {
223  auto event = (*mixHepMC_itr).GetEvent();
224 
225  size_t num_particles = event->particles_size();
226 
227  // Fill output collection
228  unsigned int number_of_protons = 0;
229  for (auto p = event->particles_begin(); p != event->particles_end(); ++p) {
230  HepMC::GenParticle const* part = *p;
231  if (select_(part, minPz_)) {
232  reco::GenParticle cand;
233  convertParticle_(cand, part);
234  ++number_of_protons;
235  cands.push_back(cand);
236  }
237  }
238  LogDebug("GenPUProtonProducer") << "Idx : " << idx_mix << " Bunch : " << bunch
239  << " Number of particles : " << num_particles
240  << " Number of protons : " << number_of_protons << endl;
241 
242  total_number_of_protons += number_of_protons;
243  }
244  }
245  LogDebug("GenPUProtonProducer") << "Total number of protons : " << total_number_of_protons << endl;
246  LogDebug("GenPUProtonProducer") << "Output collection size : " << cands.size() << endl;
247 
248  evt.put(std::move(candsPtr));
249 }
250 
252 
#define LogDebug(id)
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
std::shared_ptr< ConvertParticle > globalBeginRun(const edm::Run &, const edm::EventSetup &) const override
std::vector< int > bunchList_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > mixToken_
void setCollisionId(int s)
Definition: GenParticle.h:38
HepPDT::ParticleDataTable ParticleDataTable
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void produce(edm::StreamID, edm::Event &e, const edm::EventSetup &) const override
int init
Definition: HydjetWrapper.h:67
selection
main part
Definition: corrVsCorr.py:100
Run const & getRun() const
Definition: Event.cc:99
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void setVertex(const Point &vertex) override
set vertex
bool getData(T &iHolder) const
Definition: EventSetup.h:111
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double p4[4]
Definition: TauolaWrapper.h:92
T Abs(T a)
Definition: MathUtil.h:49
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GenPUProtonProducer(const edm::ParameterSet &)
RunIndex index() const
Definition: Run.cc:21
double f[11][100]
static int PDGCacheMax
T const * product() const
Definition: Handle.h:74
void setThreeCharge(Charge qx3) final
set electric charge
Definition: LeafCandidate.h:97
static const double mmToCm
part
Definition: HCALResponse.h:20
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
fixed size matrix
HLT enums.
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
void setStatus(int status) final
set status word
void setPdgId(int pdgId) final
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511
#define constexpr
Definition: event.py:1
Definition: Run.h:45