CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
31 
35 
36 // pythia
37 #include <Pythia8/Pythia.h>
38 
39 // root
40 #include "TH1D.h"
41 #include "TFile.h"
42 #include "TLorentzVector.h"
43 
44 //
45 // class declaration
46 //
47 
49 public:
50  explicit TestPythiaDecays(const edm::ParameterSet&);
52 
53  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
54 
55 
56 private:
57  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
58 
59  //virtual void beginRun(edm::Run const&, edm::EventSetup const&) override;
60  //virtual void endRun(edm::Run const&, edm::EventSetup const&) override;
61  //virtual void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
62  //virtual void endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
63 
64  // ----------member data ---------------------------
65  std::vector<int> pids;
66  std::map<int,TH1D*> h_mass;
67  std::map<int,TH1D*> h_p;
68  std::map<int,TH1D*> h_v;
69  std::map<int,TH1D*> h_mass_ref;
70  std::map<int,TH1D*> h_plt; // plt: proper life time
71  std::map<int,TH1D*> h_originVertexRho;
72  std::map<int,TH1D*> h_originVertexZ;
73  std::map<int,TH1D*> h_decayVertexRho;
74  std::map<int,TH1D*> h_decayVertexZ;
75  std::map<int,TH1D*> h_plt_ref; // plt: proper life time
76  std::map<int,TH1D*> h_br;
77  std::map<int,TH1D*> h_br_ref;
78 
79  std::map<int,std::vector<string> > knownDecayModes;
80 
81  Pythia8::Pythia * pythia;
83 };
84 
85 //
86 // constants, enums and typedefs
87 //
88 
89 //
90 // static data member definitions
91 //
92 
93 //
94 // constructors and destructor
95 //
97 {
98 
99  // output file
100  outputFile = iConfig.getParameter<std::string>("outputFile");
101 
102  // create pythia8 instance to access particle data
103  pythia = new Pythia8::Pythia();
104  pythia->init();
105  Pythia8::ParticleData pdt = pythia->particleData;
106 
107  // which particles will we study?
108  pids.push_back(15); // tau
109  pids.push_back(211); // pi+
110  pids.push_back(111); // pi0
111  pids.push_back(130); // K0L
112  pids.push_back(321); // K+
113  pids.push_back(323); // K*(392)
114  pids.push_back(411); // D+
115  pids.push_back(521); // B+
116 
117  // define histograms
118  for(size_t i = 0;i<pids.size();++i){
119 
120  int pid = abs(pids[i]);
121 
122  // get particle data
123  if(!pdt.isParticle(pid)){
124  std::cout << "ERROR: BAD PARTICLE, pythia is not aware of pid " << pid << std::endl;
125  std::exit(1);
126  }
127  Pythia8::ParticleDataEntry * pd = pdt.particleDataEntryPtr(pid);
128 
129  // mass histograms
130  double m0 = pd->m0();
131  double w = pd->mWidth();
132  double mmin,mmax;
133  if( w == 0){
134  mmin = m0 - m0/1000.;
135  mmax = m0 + m0/1000.;
136  }
137  else{
138  mmin = m0 - 2*w;
139  mmax = m0 + 2*w;
140  }
141  std::stringstream strstr;
142  strstr << "mass_" << pid;
143  h_mass[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,mmin,mmax);
144  h_mass_ref[pid] = (TH1D*)(h_mass[pid]->Clone(strstr.str().c_str()));
145  h_mass_ref[pid]->SetTitle(h_mass_ref[pid]->GetName());
146  if(w==0)
147  h_mass_ref[pid]->Fill(m0);
148  else{
149  for(int b =1;b<=h_mass_ref[pid]->GetNbinsX();++b){
150  double _val = h_mass_ref[pid]->GetBinCenter(b);
151  h_mass_ref[pid]->SetBinContent(b,TMath::BreitWigner(_val,m0,w));
152  }
153  }
154 
155  // p histogram
156  strstr.str("");
157  strstr << "p_" << pid;
158  h_p[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,20);
159 
160  // v histogram
161  strstr.str("");
162  strstr << "v_" << pid;
163  h_v[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,1.);
164 
165  // ctau histograms
166  double ctau0 = pd->tau0()/10.;
167  strstr.str("");
168  strstr << "plt_" << pid;
169  h_plt[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,std::min(5.*ctau0,500.));
170  h_plt_ref[pid] = (TH1D*)(h_plt[pid]->Clone(strstr.str().c_str()));
171  h_plt_ref[pid]->SetTitle(h_plt_ref[pid]->GetName());
172  for(int b =1;b<=h_plt_ref[pid]->GetNbinsX();++b){
173  double _val = h_plt_ref[pid]->GetBinCenter(b);
174  h_plt_ref[pid]->SetBinContent(b,TMath::Exp(-_val/ctau0)); //convert mm to cm
175  }
176 
177 
178  // br histograms
179  strstr.str("");
180  strstr << "br_" << pid;
181  h_br[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),0,0,0);
182  h_br[pid]->SetBit(TH1::kCanRebin);
183  h_br_ref[pid] = (TH1D*)(h_br[pid]->Clone(strstr.str().c_str()));
184  h_br_ref[pid]->SetTitle(h_br_ref[pid]->GetName());
185  knownDecayModes[pid] = vector<string>();
186  for(int d = 0;d<pd->sizeChannels();++d){
187  Pythia8::DecayChannel & channel = pd->channel(d);
188  std::vector<int> prod;
189  for(int p = 0;p<channel.multiplicity();++p){
190  int pId = abs(channel.product(p));
191  // from FastSimulation/Event/src/KineParticleFilter.cc
192  bool particleCut = ( pId > 10 && pId != 12 && pId != 14 &&
193  pId != 16 && pId != 18 && pId != 21 &&
194  (pId < 23 || pId > 40 ) &&
195  (pId < 81 || pId > 100 ) && pId != 2101 &&
196  pId != 3101 && pId != 3201 && pId != 1103 &&
197  pId != 2103 && pId != 2203 && pId != 3103 &&
198  pId != 3203 && pId != 3303 );
199  if(particleCut)
200  prod.push_back(abs(channel.product(p)));
201  }
202  std::sort(prod.begin(),prod.end());
203  strstr.str("");
204  for(size_t p = 0;p<prod.size();++p){
205  strstr << "_" << prod[p];
206  }
207  std::string label = strstr.str();
208  h_br[pid]->Fill(label.c_str(),0.);
209  h_br_ref[pid]->Fill(label.c_str(),channel.bRatio());
210  h_br[pid]->SetEntries(0);
211  knownDecayModes[pid].push_back(label);
212  }
213 
214  // vertex plots
215  strstr.str("");
216  strstr << "originVertexRho_" << pid;
217  h_originVertexRho[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,200);
218  strstr.str("");
219  strstr << "originVertexZ_" << pid;
220  h_originVertexZ[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,400);
221  strstr.str("");
222  strstr << "decayVertexRho_" << pid;
223  h_decayVertexRho[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,200);
224  strstr.str("");
225  strstr << "decayVertexZ_" << pid;
226  h_decayVertexZ[pid] = new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,400);
227  }
228 
229 }
230 
231 
233 {
234 
235  // do anything here that needs to be done at desctruction time
236  // (e.g. close files, deallocate resources etc.)
237  TFile * f = TFile::Open(outputFile.c_str(),"RECREATE");
238  f->cd();
239  f->mkdir("observed");
240  f->mkdir("prediction");
241  for(size_t i = 0;i<pids.size();++i){
242  int pid = pids[i];
243  f->cd("observed");
244  h_mass[pid]->Write();
245  h_plt[pid]->Write();
246  h_br[pid]->Write();
247  h_originVertexZ[pid]->Write();
248  h_originVertexRho[pid]->Write();
249  h_decayVertexZ[pid]->Write();
250  h_decayVertexRho[pid]->Write();
251  h_p[pid]->Write();
252  h_v[pid]->Write();
253  f->cd("prediction");
254  h_mass_ref[pid]->Write();
255  h_plt_ref[pid]->Write();
256  h_br_ref[pid]->Write();
257  }
258  f->Close();
259 }
260 
261 
262 //
263 // member functions
264 //
265 
266 // ------------ method called for each event ------------
267 void
269 {
270  using namespace edm;
271 
272  Handle<std::vector<SimTrack> > simtracks;
273  iEvent.getByLabel("famosSimHits",simtracks);
274 
275  Handle<std::vector<SimVertex> > simvertices;
276  iEvent.getByLabel("famosSimHits",simvertices);
277 
278  // create maps
279 
280  // initialize
281  std::map<size_t,std::vector<size_t> > childMap; // child indices vs parent index
282  std::map<size_t,int> parentMap; // parent index vs child index
283  for(size_t j = 0;j<simtracks->size();j++){
284  childMap[j] = std::vector<size_t>();
285  parentMap[j] = -1;
286  }
287 
288  // do the mapping
289  for(size_t j = 0;j<simtracks->size();j++){
290  size_t childIndex = j;
291  const SimTrack & child = simtracks->at(childIndex);
292  if(child.noVertex())
293  continue;
294  const SimVertex & vertex = simvertices->at(child.vertIndex());
295  if(vertex.noParent())
296  continue;
297  size_t parentIndex = vertex.parentIndex();
298  childMap[parentIndex].push_back(childIndex);
299  parentMap[childIndex] = int(parentIndex);
300  }
301 
302 
303  for(size_t j = 0;j<simtracks->size();j++){
304  const SimTrack & parent = simtracks->at(j);
305  int pid = abs(parent.type());
306  if(std::find(pids.begin(),pids.end(),pid)==pids.end())
307  continue;
308 
309  // fill mass hist
310  double mass = parent.momentum().M();
311  h_mass[pid]->Fill(mass);
312 
313  // fill p hist
314  h_p[pid]->Fill(parent.momentum().P());
315 
316  // fill vertex position hist
317  if(!parent.noVertex()){
318  const SimVertex & originVertex = simvertices->at(parent.vertIndex());
319  h_originVertexRho[pid]->Fill(originVertex.position().Rho());
320  h_originVertexZ[pid]->Fill(std::fabs(originVertex.position().Z()));
321  }
322  if(childMap[j].size() > 0){
323  const SimTrack & child = simtracks->at(childMap[j][0]);
324  const SimVertex & decayVertex = simvertices->at(child.vertIndex());
325  h_decayVertexRho[pid]->Fill(decayVertex.position().Rho());
326  h_decayVertexZ[pid]->Fill(std::fabs(decayVertex.position().Z()));
327  }
328  }
329 
330 
331  for(std::map<size_t,std::vector<size_t> >::iterator it = childMap.begin();it != childMap.end();++it){
332 
333  // fill ctau hist
334  size_t parentIndex = it->first;
335  const SimTrack & parent = simtracks->at(parentIndex);
336  int pid = abs(parent.type());
337  vector<size_t> & childIndices = it->second;
338  if(childIndices.size() == 0)
339  continue;
340 
341  if(std::find(pids.begin(),pids.end(),pid)==pids.end())
342  continue;
343 
344 
345  const SimVertex & origin_vertex = simvertices->at(parent.vertIndex());
346  const SimTrack & child0 = simtracks->at(childIndices[0]);
347  const SimVertex & decay_vertex = simvertices->at(child0.vertIndex());
348 
349  TLorentzVector lv_origin_vertex(origin_vertex.position().X(),origin_vertex.position().Y(),origin_vertex.position().Z(),origin_vertex.position().T());
350  TLorentzVector lv_decay_vertex(decay_vertex.position().X(),decay_vertex.position().Y(),decay_vertex.position().Z(),decay_vertex.position().T());
351  TLorentzVector lv_dist = lv_decay_vertex - lv_origin_vertex;
352  TLorentzVector lv_parent(parent.momentum().Px(),parent.momentum().Py(),parent.momentum().Pz(),parent.momentum().E());
353  TVector3 boost = lv_parent.BoostVector();
354  lv_dist.Boost(-boost);
355  h_v[pid]->Fill(boost.Mag());
356  double plt = lv_dist.T();
357  h_plt[pid]->Fill(plt);
358 
359  // fill br hist
360  std::vector<int> prod;
361  for(size_t d = 0;d<childIndices.size();++d){
362  prod.push_back(abs(simtracks->at(childIndices[d]).type()));
363  }
364  std::sort(prod.begin(),prod.end());
365  std::stringstream strstr;
366  for(size_t p = 0;p<prod.size();++p){
367  strstr << "_" << prod[p];
368  }
369  std::string label = strstr.str();
370  if(std::find(knownDecayModes[pid].begin(),knownDecayModes[pid].end(),label)==knownDecayModes[pid].end())
371  label = "u" + label;
372  h_br[pid]->Fill(label.c_str(),1.);
373  h_br_ref[pid]->Fill(label.c_str(),0.); // keep h_br and h_br_ref in sync
374  }
375 }
376 
377 
378 
379 // ------------ method called when starting to processes a run ------------
380 /*
381 void
382 TestPythiaDecays::beginRun(edm::Run const&, edm::EventSetup const&)
383 {
384 }
385 */
386 
387 // ------------ method called when ending the processing of a run ------------
388 /*
389 void
390 TestPythiaDecays::endRun(edm::Run const&, edm::EventSetup const&)
391 {
392 }
393 */
394 
395 // ------------ method called when starting to processes a luminosity block ------------
396 /*
397 void
398 TestPythiaDecays::beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
399 {
400 }
401 */
402 
403 // ------------ method called when ending the processing of a luminosity block ------------
404 /*
405 void
406 TestPythiaDecays::endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
407 {
408 }
409 */
410 
411 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
412 void
414  //The following says we do not know what parameters are allowed so do no validation
415  // Please change this to state exactly what you do use, even if it is no parameters
417  desc.setUnknown();
418  descriptions.addDefault(desc);
419 }
420 
421 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::map< int, std::vector< string > > knownDecayModes
list parent
Definition: dbtoconf.py:74
const double w
Definition: UKUtility.cc:23
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::map< int, TH1D * > h_v
std::vector< int > pids
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
std::map< int, TH1D * > h_decayVertexZ
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
int parentIndex() const
Definition: SimVertex.h:33
Pythia8::Pythia * pythia
std::map< int, TH1D * > h_br_ref
std::map< int, TH1D * > h_plt_ref
std::map< int, TH1D * > h_decayVertexRho
std::map< int, TH1D * > h_originVertexZ
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
bool noVertex() const
Definition: SimTrack.h:30
double f[11][100]
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:26
#define end
Definition: vmac.h:37
T min(T a, T b)
Definition: MathUtil.h:58
HepPDT::ParticleData ParticleData
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:29
std::map< int, TH1D * > h_plt
tuple pid
Definition: sysUtil.py:22
double b
Definition: hdecay.h:120
TestPythiaDecays(const edm::ParameterSet &)
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
std::string outputFile
#define begin
Definition: vmac.h:30
std::map< int, TH1D * > h_br
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
std::map< int, TH1D * > h_mass
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
tuple cout
Definition: gather_cfg.py:121
std::map< int, TH1D * > h_p
bool noParent() const
Definition: SimVertex.h:34
std::map< int, TH1D * > h_mass_ref
std::map< int, TH1D * > h_originVertexRho
tuple size
Write out results.