StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main07.cc
1 // main07.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2012 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Illustration how to generate various two-body channels from
7 // astroparticle processes, e.g. neutralino annihilation or decay.
8 // To this end a "blob" of energy is created with unit cross section,
9 // from the fictitious collision of two non-radiating incoming e+e-.
10 // In the accompanying main29.cmnd file the decay channels of this
11 // blob can be set up. Furthermore, only gamma, e+-, p/pbar and
12 // neutrinos are stable, everything else is set to decay.
13 
14 #include "Pythia.h"
15 
16 using namespace Pythia8;
17 
18 //==========================================================================
19 
20 // A derived class for (e+ e- ->) GenericResonance -> various final states.
21 
22 class Sigma1GenRes : public Sigma1Process {
23 
24 public:
25 
26  // Constructor.
27  Sigma1GenRes() {}
28 
29  // Evaluate sigmaHat(sHat): dummy unit cross section.
30  virtual double sigmaHat() {return 1.;}
31 
32  // Select flavour. No colour or anticolour.
33  virtual void setIdColAcol() {setId( -11, 11, 999999);
34  setColAcol( 0, 0, 0, 0, 0, 0);}
35 
36  // Info on the subprocess.
37  virtual string name() const {return "GenericResonance";}
38  virtual int code() const {return 9001;}
39  virtual string inFlux() const {return "ffbarSame";}
40 
41 };
42 
43 //==========================================================================
44 
45 int main() {
46 
47  // Pythia generator.
48  Pythia pythia;
49 
50  // A class to generate the fictitious resonance initial state.
51  SigmaProcess* sigma1GenRes = new Sigma1GenRes();
52 
53  // Hand pointer to Pythia.
54  pythia.setSigmaPtr( sigma1GenRes);
55 
56  // Read in the rest of the settings and data from a separate file.
57  pythia.readFile("main07.cmnd");
58 
59  // Initialization.
60  pythia.init();
61 
62  // Extract settings to be used in the main program.
63  int nEvent = pythia.mode("Main:numberOfEvents");
64  int nAbort = pythia.mode("Main:timesAllowErrors");
65 
66  // Histogram particle spectra.
67  Hist eGamma("energy spectrum of photons", 100, 0., 250.);
68  Hist eE( "energy spectrum of e+ and e-", 100, 0., 250.);
69  Hist eP( "energy spectrum of p and pbar", 100, 0., 250.);
70  Hist eNu( "energy spectrum of neutrinos", 100, 0., 250.);
71  Hist eRest( "energy spectrum of rest particles", 100, 0., 250.);
72 
73  // Begin event loop.
74  int iAbort = 0;
75  for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
76 
77  // Generate events. Quit if many failures.
78  if (!pythia.next()) {
79  if (++iAbort < nAbort) continue;
80  cout << " Event generation aborted prematurely, owing to error!\n";
81  break;
82  }
83 
84  // Loop over all particles and analyze the final-state ones.
85  for (int i = 0; i < pythia.event.size(); ++i)
86  if (pythia.event[i].isFinal()) {
87  int idAbs = pythia.event[i].idAbs();
88  double eI = pythia.event[i].e();
89  if (idAbs == 22) eGamma.fill(eI);
90  else if (idAbs == 11) eE.fill(eI);
91  else if (idAbs == 2212) eP.fill(eI);
92  else if (idAbs == 12 || idAbs == 14 || idAbs == 16) eNu.fill(eI);
93  else {
94  eRest.fill(eI);
95  cout << " Error: stable id = " << pythia.event[i].id() << endl;
96  }
97  }
98 
99  // End of event loop.
100  }
101 
102  // Final statistics and histograms.
103  pythia.stat();
104  cout << eGamma << eE << eP << eNu << eRest;
105 
106  // Done.
107  return 0;
108 }