StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main84.cc
1 // main84.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2014 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 // This program is written by Stefan Prestel.
7 // It illustrates how to do CKKW-L merging,
8 // see the Matrix Element Merging page in the online manual.
9 
10 #include <time.h>
11 #include "Pythia8/Pythia.h"
12 #include "Pythia8/Pythia8ToHepMC.h"
13 #include "HepMC/GenEvent.h"
14 #include "HepMC/IO_GenEvent.h"
15 
16 using namespace Pythia8;
17 
18 // Functions for histogramming
19 #include "fastjet/PseudoJet.hh"
20 #include "fastjet/ClusterSequence.hh"
21 #include "fastjet/CDFMidPointPlugin.hh"
22 #include "fastjet/CDFJetCluPlugin.hh"
23 #include "fastjet/D0RunIIConePlugin.hh"
24 
25 //==========================================================================
26 
27 // Find the Durham kT separation of the clustering from
28 // nJetMin --> nJetMin-1 jets in te input event
29 
30 double pTfirstJet( const Event& event, int nJetMin, double Rparam) {
31 
32  double yPartonMax = 4.;
33 
34  // Fastjet analysis - select algorithm and parameters
35  fastjet::Strategy strategy = fastjet::Best;
36  fastjet::RecombinationScheme recombScheme = fastjet::E_scheme;
37  fastjet::JetDefinition *jetDef = NULL;
38  // For hadronic collision, use hadronic Durham kT measure
39  if(event[3].colType() != 0 || event[4].colType() != 0)
40  jetDef = new fastjet::JetDefinition(fastjet::kt_algorithm, Rparam,
41  recombScheme, strategy);
42  // For e+e- collision, use e+e- Durham kT measure
43  else
44  jetDef = new fastjet::JetDefinition(fastjet::ee_kt_algorithm,
45  recombScheme, strategy);
46  // Fastjet input
47  std::vector <fastjet::PseudoJet> fjInputs;
48  // Reset Fastjet input
49  fjInputs.resize(0);
50 
51  // Loop over event record to decide what to pass to FastJet
52  for (int i = 0; i < event.size(); ++i) {
53  // (Final state && coloured+photons) only!
54  if ( !event[i].isFinal()
55  || event[i].isLepton()
56  || event[i].id() == 23
57  || abs(event[i].id()) == 24
58  || abs(event[i].y()) > yPartonMax)
59  continue;
60 
61  // Store as input to Fastjet
62  fjInputs.push_back( fastjet::PseudoJet (event[i].px(),
63  event[i].py(), event[i].pz(),event[i].e() ) );
64  }
65 
66  // Do nothing for empty input
67  if (int(fjInputs.size()) == 0) {
68  delete jetDef;
69  return 0.0;
70  }
71 
72  // Run Fastjet algorithm
73  fastjet::ClusterSequence clustSeq(fjInputs, *jetDef);
74  // Extract kT of first clustering
75  double pTFirst = sqrt(clustSeq.exclusive_dmerge_max(nJetMin-1));
76 
77  delete jetDef;
78  // Return kT
79  return pTFirst;
80 
81 }
82 
83 //==========================================================================
84 
85 int main( int argc, char* argv[] ){
86 
87  // Check that correct number of command-line arguments
88  if (argc != 6) {
89  cerr << " Unexpected number of command-line arguments. \n You are"
90  << " expected to provide the arguments" << endl
91  << " 1. Input file for settings" << endl
92  << " 2. Name of output HepMC file" << endl
93  << " 3. Maximal number of additional jets"
94  << " (not used internally in Pythia, only used to construct the full"
95  << " name of lhe files with additional jets, and to label output"
96  << " histograms)" << endl
97  << " 4. Full name of the input LHE file (with path"
98  << " , without any _0.lhe suffix)" << endl
99  << " 5. Path for output histogram files" << endl
100  << " Program stopped. " << endl;
101  return 1;
102  }
103 
104  Pythia pythia;
105 
106  // First argument: Get input from an input file
107  pythia.readFile(argv[1]);
108  int nEvent = pythia.mode("Main:numberOfEvents");
109 
110  // Interface for conversion from Pythia8::Event to HepMC event.
111  // Will fill cross section and event weight directly in this program,
112  // so switch it off for normal conversion routine.
113  HepMC::Pythia8ToHepMC ToHepMC;
114  ToHepMC.set_store_xsec(false);
115 
116  // Specify file where HepMC events will be stored.
117  HepMC::IO_GenEvent ascii_io(argv[2], std::ios::out);
118 
119  // Third argument: Maximal number of additional jets
120  int njet = atoi(argv[3]);
121 
122  // Read input and output paths
123  string iPath = string(argv[4]);
124  string oPath = string(argv[5]);
125 
126  // To write correctly normalized events to hepmc file, first get
127  // a reasonable accurate of the cross section
128  int njetCounterEstimate = njet;
129  vector<double> xsecEstimate;
130 
131  vector<double> nTrialEstimate;
132  vector<double> nAcceptEstimate;
133 
134  pythia.readString("Random:setSeed = on");
135  pythia.readString("Random:seed = 42390964");
136 
137  while(njetCounterEstimate >= 0) {
138 
139  // Number of runs
140  int nRun = 1;
141  double nTrial = 0.;
142  double nAccept = 0.;
143 
144  int countEvents = 0;
145 
146  // Run pythia nRun times with the same lhe file to get nRun times
147  // higher statistics in the histograms
148  for(int n = 1; n <= nRun ; ++n ) {
149 
150  // Get process and events from LHE file, initialize only the
151  // first time
152  bool skipInit = false;
153  if(n > 1){
154  skipInit = true;
155  pythia.readString("Main:LHEFskipInit = on");
156  }
157 
158  // From njet, choose LHE file
159  stringstream in;
160  in << "_" << njetCounterEstimate << ".lhe";
161 
162  string LHEfile = iPath + in.str();
163 
164  pythia.readString("HadronLevel:all = off");
165 
166  // Read in ME configurations
167  pythia.init(LHEfile,skipInit);
168 
169  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
170  countEvents++;
171 
172  nTrial += 1.;
173  if(iEvent == 0) pythia.stat();
174 
175  // Generate next event
176  if(pythia.next()) nAccept += 1.;
177 
178  if(countEvents == nEvent*nRun-1){
179  xsecEstimate.push_back(pythia.info.sigmaGen());
180  nTrialEstimate.push_back(nTrial+1.);
181  nAcceptEstimate.push_back(nAccept+1.);
182  }
183 
184 
185  } // end loop over events to generate
186 
187  } // end outer loop to rerun pythia with the same lhe file
188 
189  // Restart with ME of a reduced the number of jets
190  if( njetCounterEstimate > 0 )
191  njetCounterEstimate--;
192  else
193  break;
194 
195  } // end loop over different jet multiplicities
196 
197  cout << endl << "Finished estimating cross section"
198  << endl;
199 
200  for(int i=0; i < int(xsecEstimate.size()); ++i)
201  cout << " Cross section estimate for " << njet-i << " jets :"
202  << scientific << setprecision(8) << xsecEstimate[i]
203  << endl;
204  for(int i=0; i < int(nTrialEstimate.size()); ++i)
205  cout << " Trial events for " << njet-i << " jets :"
206  << scientific << setprecision(3) << nTrialEstimate[i]
207  << " Accepted events for " << njet-i << " jets :"
208  << scientific << setprecision(3) << nAcceptEstimate[i]
209  << endl;
210  cout << endl;
211 
212  // Now start merging procedure
213  int njetCounter = njet;
214 
215  Hist histPTFirstSum("pT of first jet",100,0.,100.);
216  Hist histPTSecondSum("pT of second jet",100,0.,100.);
217 
218  pythia.readString("Random:setSeed = on");
219  pythia.readString("Random:seed = 42390964");
220 
221  // Sum of event weights
222  double sigma = 0.0;
223  double sigma2 = 0.0;
224 
225  while(njetCounter >= 0) {
226 
227  cout << " Path to lhe files: " << iPath << "_*" << endl;
228  cout << " Output written to: " << oPath << "'name'.dat" << endl;
229 
230  // Set up histograms of pT of the first jet
231  Hist histPTFirst("pT of first jet",100,0.,200.);
232  Hist histPTSecond("pT of second jet",100,0.,200.);
233  Hist histPTThird("pT of third jet",100,0.,200.);
234  Hist histPTFourth("pT of fourth jet",50,0.,100.);
235  Hist histPTFifth("pT of fifth jet",30,0.,50.);
236  Hist histPTSixth("pT of sixth jet",30,0.,50.);
237 
238  // Number of runs
239  int nRun = 1;
240  // Number of tried events
241  int nTriedEvents = 0;
242  // Number of accepted events
243  int nAccepEvents = 0;
244 
245  // Run pythia nRun times with the same lhe file to get nRun times
246  // higher statistics in the histograms
247  for(int n = 1; n <= nRun ; ++n ) {
248 
249  // Get process and events from LHE file, initialize only the
250  // first time
251  bool skipInit = false;
252  if(n > 1){
253  skipInit = true;
254  pythia.readString("Main:LHEFskipInit = on");
255  }
256 
257  // From njet, choose LHE file
258  stringstream in;
259  in << "_" << njetCounter << ".lhe";
260 
261  string LHEfile = iPath + in.str();
262 
263  cout << endl << endl
264  << "\t LHE FILE FOR + " << njetCounter
265  << " JET SAMPLE READ FROM " << LHEfile
266  << endl << endl;
267 
268  cout << "Normalise with xsection " << xsecEstimate[njet-njetCounter]
269  << endl << endl;
270 
271  pythia.readString("HadronLevel:all = on");
272 
273  // Read in ME configurations
274  pythia.init(LHEfile,skipInit);
275 
276  if(pythia.flag("Main:showChangedSettings")) {
277  pythia.settings.listChanged();
278  }
279 
280  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
281 
282  nTriedEvents++;
283  if(iEvent == 0) pythia.stat();
284 
285  // Generate next event
286  if( pythia.next()) {
287 
288  double weight = pythia.info.mergingWeight();
289  nAccepEvents++;
290 
291  // Jet pT's
292  double D = 0.4;
293  double pTfirst = pTfirstJet(pythia.event,1, D);
294  double pTsecnd = pTfirstJet(pythia.event,2, D);
295  double pTthird = pTfirstJet(pythia.event,3, D);
296  double pTfourt = pTfirstJet(pythia.event,4, D);
297  double pTfifth = pTfirstJet(pythia.event,5, D);
298  double pTsixth = pTfirstJet(pythia.event,6, D);
299  histPTFirst.fill( pTfirst, weight);
300  histPTSecond.fill( pTsecnd, weight);
301  histPTThird.fill( pTthird, weight);
302  histPTFourth.fill( pTfourt, weight);
303  histPTFifth.fill( pTfifth, weight);
304  histPTSixth.fill( pTsixth, weight);
305 
306  if(weight > 0.){
307  // Construct new empty HepMC event and fill it.
308  // Units will be as chosen for HepMC build, but can be changed
309  // by arguments, e.g. GenEvt( HepMC::Units::GEV, HepMC::Units::MM)
310  HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
311 
312  double normhepmc = 1.* xsecEstimate[njet-njetCounter]
313  * nTrialEstimate[njet-njetCounter]
314  / nAcceptEstimate[njet-njetCounter]
315  * 1./ (double(nRun)*double(nEvent));
316 
317  sigma += weight*normhepmc;
318  sigma2 += pow(weight*normhepmc,2);
319  // Set event weight
320  hepmcevt->weights().push_back(weight*normhepmc);
321 
322  // Fill summed histograms
323  histPTFirstSum.fill( pTfirst, weight*normhepmc);
324  histPTSecondSum.fill( pTsecnd, weight*normhepmc);
325 
326  // Fill HepMC event, with PDF info.
327  ToHepMC.fill_next_event( pythia, hepmcevt );
328 
329  // Report cross section to hepmc
331  xsec.set_cross_section( sigma*1e9,
332  pythia.info.sigmaErr()*1e9 );
333  hepmcevt->set_cross_section( xsec );
334 
335  // Write the HepMC event to file. Done with it.
336  ascii_io << hepmcevt;
337  delete hepmcevt;
338  }
339 
340  } // if( pythia.next() )
341 
342  if(nTriedEvents%10000 == 0)
343  cout << nTriedEvents << endl;
344 
345  } // end loop over events to generate
346 
347  // print cross section, errors
348  pythia.stat();
349 
350  } // end outer loop to rerun pythia with the same lhe file
351 
352  // Normalise histograms for this particular multiplicity
353  double norm = 1.
354  * pythia.info.sigmaGen()
355  * double(nTriedEvents)/double(nAccepEvents)
356  * 1./ (double(nRun)*double(nEvent));
357 
358  histPTFirst *= norm;
359  histPTSecond *= norm;
360  histPTThird *= norm;
361  histPTFourth *= norm;
362  histPTFifth *= norm;
363  histPTSixth *= norm;
364 
365  // Write histograms for this particular multiplicity to file
366  ofstream write;
367  stringstream suffix;
368  suffix << njet << "_" << njetCounter;
369  suffix << "_wv.dat";
370 
371  write.open( (char*)(oPath + "PTjet1_" + suffix.str()).c_str());
372  histPTFirst.table(write);
373  write.close();
374 
375  write.open( (char*)(oPath + "PTjet2_" + suffix.str()).c_str());
376  histPTSecond.table(write);
377  write.close();
378 
379  write.open( (char*)(oPath + "PTjet3_" + suffix.str()).c_str());
380  histPTThird.table(write);
381  write.close();
382 
383  write.open( (char*)(oPath + "PTjet4_" + suffix.str()).c_str());
384  histPTFourth.table(write);
385  write.close();
386 
387  write.open( (char*)(oPath + "PTjet5_" + suffix.str()).c_str());
388  histPTFifth.table(write);
389  write.close();
390 
391  write.open( (char*)(oPath + "PTjet6_" + suffix.str()).c_str());
392  histPTSixth.table(write);
393  write.close();
394 
395  histPTFirst.null();
396  histPTSecond.null();
397  histPTThird.null();
398  histPTFourth.null();
399  histPTFifth.null();
400  histPTSixth.null();
401 
402  // Restart with ME of a reduced the number of jets
403  if( njetCounter > 0 )
404  njetCounter--;
405  else
406  break;
407 
408  } // end loop over different jet multiplicities
409 
410  // Since the histograms have been filled with the correct weight for
411  // each jet multiplicity, no normalisation is needed.
412  // Write summed histograms to file.
413  ofstream writeSum;
414  stringstream suffixSum;
415  suffixSum << njet << "_wv.dat";
416 
417  writeSum.open( (char*)(oPath + "PTjet1Sum_" + suffixSum.str()).c_str());
418  histPTFirstSum.table(writeSum);
419  writeSum.close();
420 
421  writeSum.open( (char*)(oPath + "PTjet2Sum_" + suffixSum.str()).c_str());
422  histPTSecondSum.table(writeSum);
423  writeSum.close();
424 
425  for(int i=0; i < int(xsecEstimate.size()); ++i)
426  cout << " Cross section estimate for " << njet-i << " jets :"
427  << scientific << setprecision(8) << xsecEstimate[i]
428  << endl;
429  for(int i=0; i < int(nTrialEstimate.size()); ++i)
430  cout << " Trial events for " << njet-i << " jets :"
431  << scientific << setprecision(3) << nTrialEstimate[i]
432  << " Accepted events for " << njet-i << " jets :"
433  << scientific << setprecision(3) << nAcceptEstimate[i]
434  << endl;
435  cout << endl;
436 
437  cout << "Histogrammed cross section for "
438  << iPath << " with " << njet << " additional jets is "
439  << scientific << setprecision(8) << sigma
440  << " error " << sqrt(sigma2) << endl;
441 
442  return 0;
443 }
void set_cross_section(double xs, double xs_err)
Set cross section and error in pb.
The GenCrossSection class stores the generated cross section.
void set_cross_section(const GenCrossSection &)
provide a pointer to the GenCrossSection container
Definition: GenEvent.h:752
The GenEvent class is the core of HepMC.
Definition: GenEvent.h:155
IO_GenEvent also deals with HeavyIon and PdfInfo.
Definition: IO_GenEvent.h:63
WeightContainer & weights()
direct access to WeightContainer
Definition: GenEvent.h:699