CMS 3D CMS Logo

TestPythiaDecays.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Lukas/TestPythiaDecays
4 // Class: TestPythiaDecays
5 //
13 //
14 // Original Author: Lukas Vanelderen
15 // Created: Tue, 13 May 2014 09:50:05 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
25 
31 
35 
39 
40 // pythia
41 #include <Pythia8/Pythia.h>
42 
43 // root
44 #include "TH1D.h"
45 #include "TLorentzVector.h"
46 
47 #if PYTHIA_VERSION_INTEGER < 8304
48 typedef Pythia8::ParticleDataEntry* ParticleDataEntryPtr;
49 #else
51 #endif
52 //
53 // class declaration
54 //
55 
56 class TestPythiaDecays : public edm::one::EDAnalyzer<edm::one::SharedResources> {
57 public:
58  explicit TestPythiaDecays(const edm::ParameterSet&);
59  ~TestPythiaDecays() override = default;
60 
61  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
62 
63 private:
64  void analyze(const edm::Event&, const edm::EventSetup&) override;
65 
66  // ----------member data ---------------------------
69  std::vector<int> pids;
70  std::map<int, TH1D*> h_mass;
71  std::map<int, TH1D*> h_p;
72  std::map<int, TH1D*> h_v;
73  std::map<int, TH1D*> h_mass_ref;
74  std::map<int, TH1D*> h_plt; // plt: proper life time
75  std::map<int, TH1D*> h_originVertexRho;
76  std::map<int, TH1D*> h_originVertexZ;
77  std::map<int, TH1D*> h_decayVertexRho;
78  std::map<int, TH1D*> h_decayVertexZ;
79  std::map<int, TH1D*> h_plt_ref; // plt: proper life time
80  std::map<int, TH1D*> h_br;
81  std::map<int, TH1D*> h_br_ref;
82 
83  std::map<int, std::vector<std::string> > knownDecayModes;
84 
85  Pythia8::Pythia* pythia;
87 };
88 
89 //
90 // constants, enums and typedefs
91 //
92 
93 //
94 // static data member definitions
95 //
96 
97 //
98 // constructors and destructor
99 //
101  : tokTrack_(consumes<std::vector<SimTrack> >(edm::InputTag("famosSimHits"))),
102  tokVertex_(consumes<std::vector<SimVertex> >(edm::InputTag("famosSimHits"))) {
103  usesResource(TFileService::kSharedResource);
104  // output file
105  outputFile = iConfig.getParameter<std::string>("outputFile");
106 
107  // create pythia8 instance to access particle data
108  pythia = new Pythia8::Pythia();
109  pythia->init();
110  Pythia8::ParticleData pdt = pythia->particleData;
111 
112  // which particles will we study?
113  pids.push_back(15); // tau
114  pids.push_back(211); // pi+
115  pids.push_back(111); // pi0
116  pids.push_back(130); // K0L
117  pids.push_back(321); // K+
118  pids.push_back(323); // K*(392)
119  pids.push_back(411); // D+
120  pids.push_back(521); // B+
121 
122  // define histograms
124  TFileDirectory observe = file->mkdir("observed");
125  TFileDirectory predict = file->mkdir("prediction");
126  for (size_t i = 0; i < pids.size(); ++i) {
127  int pid = std::abs(pids[i]);
128 
129  // get particle data
130  if (!pdt.isParticle(pid)) {
131  edm::LogError("TestPythiaDecays") << "ERROR: BAD PARTICLE, pythia is not aware of pid " << pid;
132  throw cms::Exception("Unknown", "FastSim") << "Bad Particle Type " << pid << "\n";
133  std::exit(1);
134  }
135  ParticleDataEntryPtr pd = pdt.particleDataEntryPtr(pid);
136 
137  // mass histograms
138  double m0 = pd->m0();
139  double w = pd->mWidth();
140  double mmin, mmax;
141  if (w == 0) {
142  mmin = m0 - m0 / 1000.;
143  mmax = m0 + m0 / 1000.;
144  } else {
145  mmin = m0 - 2 * w;
146  mmax = m0 + 2 * w;
147  }
148  std::stringstream strstr;
149  strstr << "mass_" << pid;
150  h_mass[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, mmin, mmax);
151  h_mass_ref[pid] = predict.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, mmin, mmax);
152  if (w == 0)
153  h_mass_ref[pid]->Fill(m0);
154  else {
155  for (int b = 1; b <= h_mass_ref[pid]->GetNbinsX(); ++b) {
156  double _val = h_mass_ref[pid]->GetBinCenter(b);
157  h_mass_ref[pid]->SetBinContent(b, TMath::BreitWigner(_val, m0, w));
158  }
159  }
160 
161  // p histogram
162  strstr.str("");
163  strstr << "p_" << pid;
164  h_p[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, 20);
165 
166  // v histogram
167  strstr.str("");
168  strstr << "v_" << pid;
169  h_v[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, 1.);
170 
171  // ctau histograms
172  double ctau0 = pd->tau0() / 10.;
173  strstr.str("");
174  strstr << "plt_" << pid;
175  h_plt[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, std::min(5. * ctau0, 500.));
176  h_plt_ref[pid] = predict.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, std::min(5. * ctau0, 500.));
177  h_plt_ref[pid]->SetTitle(h_plt_ref[pid]->GetName());
178  for (int b = 1; b <= h_plt_ref[pid]->GetNbinsX(); ++b) {
179  double _val = h_plt_ref[pid]->GetBinCenter(b);
180  h_plt_ref[pid]->SetBinContent(b, TMath::Exp(-_val / ctau0)); //convert mm to cm
181  }
182 
183  // br histograms
184  strstr.str("");
185  strstr << "br_" << pid;
186  h_br[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 0, 0, 0);
187  h_br[pid]->SetCanExtend(TH1::kAllAxes);
188  h_br_ref[pid] = predict.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 0, 0, 0);
189  h_br_ref[pid]->SetCanExtend(TH1::kAllAxes);
190  h_br_ref[pid]->SetTitle(h_br_ref[pid]->GetName());
191  knownDecayModes[pid] = std::vector<std::string>();
192  for (int d = 0; d < pd->sizeChannels(); ++d) {
193  Pythia8::DecayChannel& channel = pd->channel(d);
194  std::vector<int> prod;
195  for (int p = 0; p < channel.multiplicity(); ++p) {
196  int pId = abs(channel.product(p));
197  // from FastSimulation/Event/src/KineParticleFilter.cc
198  bool particleCut =
199  (pId > 10 && pId != 12 && pId != 14 && pId != 16 && pId != 18 && pId != 21 && (pId < 23 || pId > 40) &&
200  (pId < 81 || pId > 100) && pId != 2101 && pId != 3101 && pId != 3201 && pId != 1103 && pId != 2103 &&
201  pId != 2203 && pId != 3103 && pId != 3203 && pId != 3303);
202  if (particleCut)
203  prod.push_back(abs(channel.product(p)));
204  }
205  std::sort(prod.begin(), prod.end());
206  strstr.str("");
207  for (size_t p = 0; p < prod.size(); ++p) {
208  strstr << "_" << prod[p];
209  }
210  std::string label = strstr.str();
211  h_br[pid]->Fill(label.c_str(), 0.);
212  h_br_ref[pid]->Fill(label.c_str(), channel.bRatio());
213  h_br[pid]->SetEntries(0);
214  knownDecayModes[pid].push_back(label);
215  }
216 
217  // vertex plots
218  strstr.str("");
219  strstr << "originVertexRho_" << pid;
220  h_originVertexRho[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, 200);
221  strstr.str("");
222  strstr << "originVertexZ_" << pid;
223  h_originVertexZ[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, 400);
224  strstr.str("");
225  strstr << "decayVertexRho_" << pid;
226  h_decayVertexRho[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, 200);
227  strstr.str("");
228  strstr << "decayVertexZ_" << pid;
229  h_decayVertexZ[pid] = observe.make<TH1D>(strstr.str().c_str(), strstr.str().c_str(), 100, 0, 400);
230  }
231 }
232 
233 //
234 // member functions
235 //
236 
237 // ------------ method called for each event ------------
239  const edm::Handle<std::vector<SimTrack> >& simtracks = iEvent.getHandle(tokTrack_);
240  const edm::Handle<std::vector<SimVertex> > simvertices = iEvent.getHandle(tokVertex_);
241 
242  // create maps
243 
244  // initialize
245  std::map<size_t, std::vector<size_t> > childMap; // child indices vs parent index
246  std::map<size_t, int> parentMap; // parent index vs child index
247  for (size_t j = 0; j < simtracks->size(); j++) {
248  childMap[j] = std::vector<size_t>();
249  parentMap[j] = -1;
250  }
251 
252  // do the mapping
253  for (size_t j = 0; j < simtracks->size(); j++) {
254  size_t childIndex = j;
255  const SimTrack& child = simtracks->at(childIndex);
256  if (child.noVertex())
257  continue;
258  const SimVertex& vertex = simvertices->at(child.vertIndex());
259  if (vertex.noParent())
260  continue;
261  size_t parentIndex = vertex.parentIndex();
262  childMap[parentIndex].push_back(childIndex);
263  parentMap[childIndex] = int(parentIndex);
264  }
265 
266  for (size_t j = 0; j < simtracks->size(); j++) {
267  const SimTrack& parent = simtracks->at(j);
268  int pid = abs(parent.type());
269  if (std::find(pids.begin(), pids.end(), pid) == pids.end())
270  continue;
271 
272  // fill mass hist
273  double mass = parent.momentum().M();
274  h_mass[pid]->Fill(mass);
275 
276  // fill p hist
277  h_p[pid]->Fill(parent.momentum().P());
278 
279  // fill vertex position hist
280  if (!parent.noVertex()) {
281  const SimVertex& originVertex = simvertices->at(parent.vertIndex());
282  h_originVertexRho[pid]->Fill(originVertex.position().Rho());
283  h_originVertexZ[pid]->Fill(std::fabs(originVertex.position().Z()));
284  }
285  if (!childMap[j].empty()) {
286  const SimTrack& child = simtracks->at(childMap[j][0]);
287  const SimVertex& decayVertex = simvertices->at(child.vertIndex());
288  h_decayVertexRho[pid]->Fill(decayVertex.position().Rho());
289  h_decayVertexZ[pid]->Fill(std::fabs(decayVertex.position().Z()));
290  }
291  }
292 
293  for (std::map<size_t, std::vector<size_t> >::iterator it = childMap.begin(); it != childMap.end(); ++it) {
294  // fill ctau hist
295  size_t parentIndex = it->first;
296  const SimTrack& parent = simtracks->at(parentIndex);
297  int pid = abs(parent.type());
298  std::vector<size_t>& childIndices = it->second;
299  if (childIndices.empty())
300  continue;
301 
302  if (std::find(pids.begin(), pids.end(), pid) == pids.end())
303  continue;
304 
305  const SimVertex& origin_vertex = simvertices->at(parent.vertIndex());
306  const SimTrack& child0 = simtracks->at(childIndices[0]);
307  const SimVertex& decay_vertex = simvertices->at(child0.vertIndex());
308 
309  TLorentzVector lv_origin_vertex(origin_vertex.position().X(),
310  origin_vertex.position().Y(),
311  origin_vertex.position().Z(),
312  origin_vertex.position().T());
313  TLorentzVector lv_decay_vertex(decay_vertex.position().X(),
314  decay_vertex.position().Y(),
315  decay_vertex.position().Z(),
316  decay_vertex.position().T());
317  TLorentzVector lv_dist = lv_decay_vertex - lv_origin_vertex;
318  TLorentzVector lv_parent(
319  parent.momentum().Px(), parent.momentum().Py(), parent.momentum().Pz(), parent.momentum().E());
320  TVector3 boost = lv_parent.BoostVector();
321  lv_dist.Boost(-boost);
322  h_v[pid]->Fill(boost.Mag());
323  double plt = lv_dist.T();
324  h_plt[pid]->Fill(plt);
325 
326  // fill br hist
327  std::vector<int> prod;
328  for (size_t d = 0; d < childIndices.size(); ++d) {
329  prod.push_back(abs(simtracks->at(childIndices[d]).type()));
330  }
331  std::sort(prod.begin(), prod.end());
332  std::stringstream strstr;
333  for (size_t p = 0; p < prod.size(); ++p) {
334  strstr << "_" << prod[p];
335  }
336  std::string label = strstr.str();
337  if (std::find(knownDecayModes[pid].begin(), knownDecayModes[pid].end(), label) == knownDecayModes[pid].end())
338  label = "u" + label;
339  h_br[pid]->Fill(label.c_str(), 1.);
340  h_br_ref[pid]->Fill(label.c_str(), 0.); // keep h_br and h_br_ref in sync
341  }
342 }
343 
344 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
346  //The following says we do not know what parameters are allowed so do no validation
347  // Please change this to state exactly what you do use, even if it is no parameters
349  desc.setUnknown();
350  descriptions.addDefault(desc);
351 }
352 
353 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
std::map< int, TH1D * > h_plt
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
Definition: CLHEP.h:16
T w() const
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:21
std::vector< int > pids
std::map< int, TH1D * > h_originVertexRho
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::map< int, TH1D * > h_plt_ref
std::map< int, TH1D * > h_br_ref
std::map< int, TH1D * > h_br
void analyze(const edm::Event &, const edm::EventSetup &) override
std::map< int, TH1D * > h_originVertexZ
const edm::EDGetTokenT< std::vector< SimVertex > > tokVertex_
T * make(const Args &...args) const
make new ROOT object
char const * label
std::map< int, TH1D * > h_p
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:33
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
Pythia8::Pythia * pythia
std::map< int, TH1D * > h_mass
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HepPDT::ParticleData ParticleData
d
Definition: ztail.py:151
std::map< int, TH1D * > h_v
std::map< int, std::vector< std::string > > knownDecayModes
std::map< int, TH1D * > h_decayVertexZ
double b
Definition: hdecay.h:120
TestPythiaDecays(const edm::ParameterSet &)
std::map< int, TH1D * > h_decayVertexRho
std::string outputFile
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::map< int, TH1D * > h_mass_ref
Pythia8::ParticleDataEntry * ParticleDataEntryPtr
~TestPythiaDecays() override=default
const edm::EDGetTokenT< std::vector< SimTrack > > tokTrack_
def exit(msg="")