This file is indexed.

/usr/share/pythia8-examples/examples/main84.cc is in pythia8-examples 8.1.80-1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
// main84.cc is a part of the PYTHIA event generator.
// Copyright (C) 2013 Torbjorn Sjostrand.
// PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
// Please respect the MCnet Guidelines, see GUIDELINES for details.

// This program is written by Stefan Prestel.
// It illustrates how to do CKKW-L merging, 
// see the Matrix Element Merging page in the online manual. 

#include <time.h>
#include "Pythia8/Pythia.h"
#include "Pythia8/Pythia8ToHepMC.h"
#include "HepMC/GenEvent.h"
#include "HepMC/IO_GenEvent.h"

using namespace Pythia8;

// Functions for histogramming
#include "fastjet/PseudoJet.hh"
#include "fastjet/ClusterSequence.hh"
#include "fastjet/CDFMidPointPlugin.hh"
#include "fastjet/CDFJetCluPlugin.hh"
#include "fastjet/D0RunIIConePlugin.hh"

//==========================================================================

// Find the Durham kT separation of the clustering from
// nJetMin --> nJetMin-1 jets in te input event  

double pTfirstJet( const Event& event, int nJetMin, double Rparam) {

  double yPartonMax = 4.;

  // Fastjet analysis - select algorithm and parameters
  fastjet::Strategy               strategy = fastjet::Best;
  fastjet::RecombinationScheme    recombScheme = fastjet::E_scheme;
  fastjet::JetDefinition         *jetDef = NULL;
  // For hadronic collision, use hadronic Durham kT measure
  if(event[3].colType() != 0 || event[4].colType() != 0)
    jetDef = new fastjet::JetDefinition(fastjet::kt_algorithm, Rparam,
                                      recombScheme, strategy);
  // For e+e- collision, use e+e- Durham kT measure
  else
    jetDef = new fastjet::JetDefinition(fastjet::ee_kt_algorithm,
                                      recombScheme, strategy);
  // Fastjet input
  std::vector <fastjet::PseudoJet> fjInputs;
  // Reset Fastjet input
  fjInputs.resize(0);

  // Loop over event record to decide what to pass to FastJet
  for (int i = 0; i < event.size(); ++i) {
    // (Final state && coloured+photons) only!
    if ( !event[i].isFinal()
      || event[i].isLepton()
      || event[i].id() == 23
      || abs(event[i].id()) == 24
      || abs(event[i].y()) > yPartonMax)
      continue;

    // Store as input to Fastjet
    fjInputs.push_back( fastjet::PseudoJet (event[i].px(),
            event[i].py(), event[i].pz(),event[i].e() ) );
  }

  // Do nothing for empty input 
  if (int(fjInputs.size()) == 0) {
    delete jetDef;
    return 0.0;
  }

  // Run Fastjet algorithm
  fastjet::ClusterSequence clustSeq(fjInputs, *jetDef);
  // Extract kT of first clustering
  double pTFirst = sqrt(clustSeq.exclusive_dmerge_max(nJetMin-1));

  delete jetDef;
  // Return kT
  return pTFirst;

}

//==========================================================================

int main( int argc, char* argv[] ){

  // Check that correct number of command-line arguments
  if (argc != 6) {
    cerr << " Unexpected number of command-line arguments. \n You are"
         << " expected to provide the arguments" << endl
         << " 1. Input file for settings" << endl
         << " 2. Name of output HepMC file" << endl
         << " 3. Maximal number of additional jets"
         << " (not used internally in Pythia, only used to construct the full"
         << " name of lhe files with additional jets, and to label output"
         << " histograms)" << endl
         << " 4. Full name of the input LHE file (with path"
         << " , without any _0.lhe suffix)" << endl
         << " 5. Path for output histogram files" << endl
         << " Program stopped. " << endl;
    return 1;
  }

  Pythia pythia;

  // First argument: Get input from an input file
  pythia.readFile(argv[1]);
  int nEvent = pythia.mode("Main:numberOfEvents");

  // Interface for conversion from Pythia8::Event to HepMC event. 
  // Will fill cross section and event weight directly in this program,
  // so switch it off for normal conversion routine.
  HepMC::Pythia8ToHepMC ToHepMC;
  ToHepMC.set_store_xsec(false);

  // Specify file where HepMC events will be stored.
  HepMC::IO_GenEvent ascii_io(argv[2], std::ios::out);

  // Third argument: Maximal number of additional jets
  int njet = atoi(argv[3]);

  // Read input and output paths
  string iPath = string(argv[4]);
  string oPath = string(argv[5]);

  // To write correctly normalized events to hepmc file, first get
  // a reasonable accurate of the cross section
  int njetCounterEstimate = njet;
  vector<double> xsecEstimate;

  vector<double> nTrialEstimate;
  vector<double> nAcceptEstimate;

  pythia.readString("Random:setSeed = on");
  pythia.readString("Random:seed = 42390964");

  while(njetCounterEstimate >= 0) {

    // Number of runs
    int nRun = 1;
    double nTrial = 0.;
    double nAccept = 0.;

    int countEvents = 0;

    // Run pythia nRun times with the same lhe file to get nRun times
    // higher statistics in the histograms
    for(int n = 1; n <= nRun ; ++n ) {

      // Get process and events from LHE file, initialize only the
      // first time
      bool skipInit = false;
      if(n > 1){
        skipInit = true;
        pythia.readString("Main:LHEFskipInit = on");
      }

      // From njet, choose LHE file
      stringstream in;
      in   << "_" << njetCounterEstimate << ".lhe";

      string LHEfile = iPath + in.str();

      pythia.readString("HadronLevel:all = off");

      // Read in ME configurations
      pythia.init(LHEfile,skipInit);

      for( int iEvent=0; iEvent<nEvent; ++iEvent ){
        countEvents++;

        nTrial += 1.;
        if(iEvent == 0) pythia.stat();

        // Generate next event
        if(pythia.next()) nAccept += 1.;

        if(countEvents == nEvent*nRun-1){
          xsecEstimate.push_back(pythia.info.sigmaGen());
          nTrialEstimate.push_back(nTrial+1.);
          nAcceptEstimate.push_back(nAccept+1.);
        }


      } // end loop over events to generate

    } // end outer loop to rerun pythia with the same lhe file

    // Restart with ME of a reduced the number of jets
    if( njetCounterEstimate > 0 )
      njetCounterEstimate--;
    else
      break;

  } // end loop over different jet multiplicities

  cout << endl << "Finished estimating cross section"
    << endl;

  for(int i=0; i < int(xsecEstimate.size()); ++i)
    cout << "  Cross section estimate for " << njet-i << " jets :"
      << scientific << setprecision(8) << xsecEstimate[i]
      << endl;
  for(int i=0; i < int(nTrialEstimate.size()); ++i)
    cout << "  Trial events for " << njet-i << " jets :"
      << scientific << setprecision(3) << nTrialEstimate[i]
      << "  Accepted events for " << njet-i << " jets :"
      << scientific << setprecision(3) << nAcceptEstimate[i]
      << endl;
  cout << endl;

  // Now start merging procedure
  int njetCounter = njet;

  Hist histPTFirstSum("pT of first jet",100,0.,100.);
  Hist histPTSecondSum("pT of second jet",100,0.,100.);

  pythia.readString("Random:setSeed = on");
  pythia.readString("Random:seed = 42390964");

  // Sum of event weights
  double sigma = 0.0;
  double sigma2 = 0.0;

  while(njetCounter >= 0) {

    cout << "   Path to lhe files: " << iPath << "_*" << endl;
    cout << "   Output written to: " << oPath << "'name'.dat" << endl;

    // Set up histograms of pT of the first jet
    Hist histPTFirst("pT of first jet",100,0.,200.);
    Hist histPTSecond("pT of second jet",100,0.,200.);
    Hist histPTThird("pT of third jet",100,0.,200.);
    Hist histPTFourth("pT of fourth jet",50,0.,100.);
    Hist histPTFifth("pT of fifth jet",30,0.,50.);
    Hist histPTSixth("pT of sixth jet",30,0.,50.);

    // Number of runs
    int nRun = 1;
    // Number of tried events
    int nTriedEvents = 0;
    // Number of accepted events
    int nAccepEvents = 0;

    // Run pythia nRun times with the same lhe file to get nRun times
    // higher statistics in the histograms
    for(int n = 1; n <= nRun ; ++n ) {

      // Get process and events from LHE file, initialize only the
      // first time
      bool skipInit = false;
      if(n > 1){
        skipInit = true;
        pythia.readString("Main:LHEFskipInit = on");
      }

      // From njet, choose LHE file
      stringstream in;
      in   << "_" << njetCounter << ".lhe";

      string LHEfile = iPath + in.str();

      cout << endl << endl
        << "\t LHE FILE FOR + " << njetCounter
        << " JET SAMPLE READ FROM " << LHEfile
        << endl << endl;

      cout << "Normalise with xsection " << xsecEstimate[njet-njetCounter]
        << endl << endl;

      pythia.readString("HadronLevel:all = on");

      // Read in ME configurations
      pythia.init(LHEfile,skipInit);

      if(pythia.flag("Main:showChangedSettings")) {
        pythia.settings.listChanged();
      }

      for( int iEvent=0; iEvent<nEvent; ++iEvent ){

        nTriedEvents++;
        if(iEvent == 0) pythia.stat();

        // Generate next event
        if( pythia.next()) {

          double weight = pythia.info.mergingWeight();
          nAccepEvents++;

          // Jet pT's
          double D = 0.4;
          double pTfirst = pTfirstJet(pythia.event,1, D);
          double pTsecnd = pTfirstJet(pythia.event,2, D);
          double pTthird = pTfirstJet(pythia.event,3, D);
          double pTfourt = pTfirstJet(pythia.event,4, D);
          double pTfifth = pTfirstJet(pythia.event,5, D);
          double pTsixth = pTfirstJet(pythia.event,6, D);
          histPTFirst.fill( pTfirst, weight);
          histPTSecond.fill( pTsecnd, weight);
          histPTThird.fill( pTthird, weight);
          histPTFourth.fill( pTfourt, weight);
          histPTFifth.fill( pTfifth, weight);
          histPTSixth.fill( pTsixth, weight);

          if(weight > 0.){
            // Construct new empty HepMC event and fill it.
            // Units will be as chosen for HepMC build, but can be changed
            // by arguments, e.g. GenEvt( HepMC::Units::GEV, HepMC::Units::MM)  
            HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();

            double normhepmc = 1.* xsecEstimate[njet-njetCounter]
                * nTrialEstimate[njet-njetCounter]
                / nAcceptEstimate[njet-njetCounter]
                * 1./ (double(nRun)*double(nEvent));

            sigma += weight*normhepmc;
            sigma2 += pow(weight*normhepmc,2);
            // Set event weight
            hepmcevt->weights().push_back(weight*normhepmc);

            // Fill summed histograms
            histPTFirstSum.fill( pTfirst, weight*normhepmc);
            histPTSecondSum.fill( pTsecnd, weight*normhepmc);

            // Fill HepMC event, with PDF info.
            ToHepMC.fill_next_event( pythia, hepmcevt );

            // Report cross section to hepmc
            HepMC::GenCrossSection xsec;
            xsec.set_cross_section( sigma*1e9,
              pythia.info.sigmaErr()*1e9 );
            hepmcevt->set_cross_section( xsec );

            // Write the HepMC event to file. Done with it.
            ascii_io << hepmcevt;
            delete hepmcevt;
          }

        } // if( pythia.next() )

        if(nTriedEvents%10000 == 0)
          cout << nTriedEvents << endl;

      } // end loop over events to generate

      // print cross section, errors
      pythia.stat();

    } // end outer loop to rerun pythia with the same lhe file

    // Normalise histograms for this particular multiplicity
    double norm = 1.
                * pythia.info.sigmaGen()
                * double(nTriedEvents)/double(nAccepEvents)
                * 1./ (double(nRun)*double(nEvent));

    histPTFirst           *= norm;
    histPTSecond          *= norm;
    histPTThird           *= norm;
    histPTFourth          *= norm;
    histPTFifth           *= norm;
    histPTSixth           *= norm;

    // Write histograms for this particular multiplicity to file
    ofstream write;
    stringstream suffix;
    suffix << njet << "_" << njetCounter;
    suffix << "_wv.dat";

    write.open( (char*)(oPath + "PTjet1_" + suffix.str()).c_str());
    histPTFirst.table(write);
    write.close();

    write.open( (char*)(oPath + "PTjet2_" + suffix.str()).c_str());
    histPTSecond.table(write);
    write.close();

    write.open( (char*)(oPath + "PTjet3_" + suffix.str()).c_str());
    histPTThird.table(write);
    write.close();

    write.open( (char*)(oPath + "PTjet4_" + suffix.str()).c_str());
    histPTFourth.table(write);
    write.close();

    write.open( (char*)(oPath + "PTjet5_" + suffix.str()).c_str());
    histPTFifth.table(write);
    write.close();

    write.open( (char*)(oPath + "PTjet6_" + suffix.str()).c_str());
    histPTSixth.table(write);
    write.close();

    histPTFirst.null();
    histPTSecond.null();
    histPTThird.null();
    histPTFourth.null();
    histPTFifth.null();
    histPTSixth.null();

    // Restart with ME of a reduced the number of jets
    if( njetCounter > 0 )
      njetCounter--;
    else
      break;

  } // end loop over different jet multiplicities

  // Since the histograms have been filled with the correct weight for
  // each jet multiplicity, no normalisation is needed.
  // Write summed histograms to file.
  ofstream writeSum;
  stringstream suffixSum;
  suffixSum << njet << "_wv.dat";

  writeSum.open( (char*)(oPath + "PTjet1Sum_" + suffixSum.str()).c_str());
  histPTFirstSum.table(writeSum);
  writeSum.close();

  writeSum.open( (char*)(oPath + "PTjet2Sum_" + suffixSum.str()).c_str());
  histPTSecondSum.table(writeSum);
  writeSum.close();

  for(int i=0; i < int(xsecEstimate.size()); ++i)
    cout << "  Cross section estimate for " << njet-i << " jets :"
      << scientific << setprecision(8) << xsecEstimate[i]
      << endl;
  for(int i=0; i < int(nTrialEstimate.size()); ++i)
    cout << "  Trial events for " << njet-i << " jets :"
      << scientific << setprecision(3) << nTrialEstimate[i]
      << "  Accepted events for " << njet-i << " jets :"
      << scientific << setprecision(3) << nAcceptEstimate[i]
      << endl;
  cout << endl;

  cout << "Histogrammed cross section for "
     << iPath << " with " << njet << " additional jets is " 
     << scientific << setprecision(8) << sigma
     << " error " << sqrt(sigma2) << endl;

  return 0;
}