StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main28.cc
1 // main28.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2012 Peter Skands, 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 // Example of of R-hadron production.
7 // Several of the possibilities shown here, like displaced vertices,
8 // are extras that need not be used for the basic setup.
9 
10 #include "Pythia.h"
11 
12 using namespace Pythia8;
13 
14 int main() {
15 
16  // Key settings to be used in the main program.
17  // nGluino = 0, 1, 2 give stop pair, single gluino or gluino pair.
18  int nGluino = 2;
19  int nEvent = 200;
20  int nAbort = 3;
21  int nList = 0;
22  double eCM = 7000.;
23 
24  // Generator. Shorthand for the event.
25  Pythia pythia;
26  Event& event = pythia.event;
27 
28  // Set up beams: p p is default so only need set energy.
29  pythia.settings.parm("Beams:eCM", eCM);
30 
31  // Squark pair: use stop-antistop as example.
32  if (nGluino == 0) {
33  pythia.readString("SUSY:gg2squarkantisquark = on");
34  pythia.readString("SUSY:idA = 1000006");
35  pythia.readString("SUSY:idB = 1000006");
36  // Squark-gluino pair: also supersymmetric u has been made long-lived.
37  // Stop does not work since then one would need inoming top PDF.
38  // Nevertheless R-hadrons are numbered/named as if containing a stop.
39  } else if (nGluino == 1) {
40  pythia.readString("SUSY:qg2squarkgluino = on");
41  pythia.readString("SUSY:idA = 1000002");
42  pythia.readString("RHadrons:idStop = 1000002");
43  pythia.readString("SUSY:idB = 1000021");
44  // Gluino pair.
45  } else {
46  pythia.readString("SUSY:gg2gluinogluino = on");
47  }
48 
49  // Use hacked sps1a file, with stop (+su) and gluino made long-lived.
50  // This is based on the width being less than 0.2 GeV by default.
51  pythia.readString("SLHA:file = sps1aNarrowStopGluino.spc");
52 
53  // Allow R-hadron formation.
54  pythia.readString("Rhadrons:allow = on");
55 
56  // If you want to do the decay separately later,
57  // you need to switch off automatic decays.
58  pythia.readString("RHadrons:allowDecay = off");
59 
60  // Fraction of gluinoballs.
61  pythia.readString("RHadrons:probGluinoball = 0.1");
62 
63  // Switch off key components.
64  //pythia.readString("PartonLevel:MPI = off");
65  //pythia.readString("PartonLevel:ISR = off");
66  //pythia.readString("PartonLevel:FSR = off");
67  //pythia.readString("HadronLevel:Hadronize = off");
68 
69  // Allow the R-hadrons to have secondary vertices: set c*tau in mm.
70  // Note that width and lifetime can be set independently.
71  // (Nonzero small widths are needed e.g. to select branching ratios.)
72  pythia.readString("1000002:tau0 = 200.");
73  pythia.readString("1000006:tau0 = 250.");
74  pythia.readString("1000021:tau0 = 300.");
75 
76  // Checks. Optionally relax E-p-conservation.
77  pythia.readString("Check:nErrList = 2");
78  //pythia.readString("Check:epTolErr = 2e-3");
79 
80  // Possibility to switch off particle data and event listings.
81  // Also to shop location of displaced vertices.
82  pythia.readString("Init:showChangedSettings = on");
83  pythia.readString("Init:showChangedParticleData = off");
84  pythia.readString("Next:numberShowInfo = 1");
85  pythia.readString("Next:numberShowProcess = 1");
86  pythia.readString("Next:numberShowEvent = 0");
87  pythia.readString("Next:showScaleAndVertex = on");
88 
89  // Initialize.
90  pythia.init();
91 
92  // Histograms.
93  Hist nChargedH("charged multiplicity", 100, -0.5, 799.5);
94  Hist dndyChargedH("dn/dy charged", 100, -10., 10.);
95  Hist dndyRH("dn/dy R-hadrons", 100, -5., 5.);
96  Hist pTRH("pT R-hadrons", 100, 0., 1000.);
97  Hist xRH("p_RHadron / p_sparticle", 100, 0.9, 1.1);
98  Hist mDiff("m(Rhadron) - m(sparticle)", 100, 0., 5.);
99  Hist decVtx("R-hadron decay vertex (mm from origin)", 100, 0., 1000.);
100 
101  // R-hadron flavour composition.
102  map<int, int> flavours;
103 
104  // Begin event loop.
105  int iAbort = 0;
106  for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
107 
108  // Generate events. Quit if failure.
109  if (!pythia.next()) {
110  if (++iAbort < nAbort) continue;
111  cout << " Event generation aborted prematurely, owing to error!\n";
112  break;
113  }
114 
115  // Loop over final charged particles in the event.
116  // The R-hadrons may not yet have decayed here.
117  int nCharged = 0;
118  Vec4 pSum;
119  for (int i = 0; i < event.size(); ++i) {
120  if (event[i].isFinal()) {
121  pSum += event[i].p();
122  if (event[i].isCharged()) {
123  ++nCharged;
124  dndyChargedH.fill( event[i].y() );
125  }
126  }
127  }
128  nChargedH.fill( nCharged );
129 
130  // Loop over final R-hadrons in the event: kinematic distribution
131  for (int i = 0; i < event.size(); ++i) {
132  int idAbs = event[i].idAbs();
133  if (idAbs > 1000100 && idAbs < 2000000 && idAbs != 1009002) {
134  ++flavours[ event[i].id() ];
135  dndyRH.fill( event[i].y() );
136  pTRH.fill( event[i].pT() );
137  // Trace back to mother; compare momenta and masses.
138  int iMother = i;
139  while( event[iMother].statusAbs() > 100)
140  iMother = event[iMother].mother1();
141  double xFrac = event[i].pAbs() / event[iMother].pAbs();
142  xRH.fill( xFrac);
143  double mShift = event[i].m() - event[iMother].m();
144  mDiff.fill( mShift );
145  // Separation of R-hadron decay vertex from origin.
146  // Don't be fooled by pAbs(); it gives the three-vector length
147  // of any Vec4, also one representing spatial coordinates.
148  double dist = event[i].vDec().pAbs();
149  decVtx.fill( dist);
150 
151  // This is a place where you could allow a R-hadron shift of
152  // identity, momentum and decay vertex to allow for detector effects.
153  // Identity not illustrated here; requires a change of mass as well.
154  // Toy model: assume an exponential energy loss, < > = 1 GeV,
155  // but at most half of kinetic energy. Unchanged direction.
156  // Note that event will no longer conserve energy and momentum.
157  double eLossAvg = 1.;
158  double eLoss = 0.;
159  do { eLoss = eLossAvg * pythia.rndm.exp(); }
160  while (eLoss > 0.5 * (event[i].e() - event[i].m()));
161  double eNew = event[i].e() - eLoss;
162  Vec4 pNew = event[i].p() * sqrt( pow2(eNew) - pow2(event[i].m()) )
163  / event[i].pAbs();
164  pNew.e( eNew);
165  event[i].p( pNew);
166  // The decay vertex will be calculated based on the production vertex,
167  // the proper lifetime tau and the NEW four-momentum, rather than
168  // e.g. some average momentum, if you do not set it by hand.
169  // This commented-out piece illustrates brute-force setting,
170  // but you should provide real numbers from some tracking program.
171  // With tau = 0 the decay is right at the chosen point.
172  //event[i].tau( 0.);
173  //event[i].vProd( 132., 155., 233., 177.);
174 
175  // End of loop over final R-hadrons.
176  }
177  }
178 
179  // If you have set R-hadrons stable above,
180  // you can still force them to decay at this stage.
181  pythia.forceRHadronDecays();
182  if (iEvent < nList) pythia.event.list(true);
183 
184  // End of event loop.
185  }
186 
187  // Final statistics, flavour composition and histogram output.
188  pythia.stat();
189  cout << "\n Composition of produced R-hadrons \n code "
190  << "name times " << endl;
191  for (map<int, int>::iterator flavNow = flavours.begin();
192  flavNow != flavours.end(); ++flavNow) cout << setw(8)
193  << flavNow->first << setw(16) << pythia.particleData.name(flavNow->first)
194  << setw(8) << flavNow->second << endl;
195  cout << nChargedH << dndyChargedH << dndyRH << pTRH << xRH << mDiff
196  << decVtx;
197 
198  // Done.
199  return 0;
200 }