File: fastjet_example.cc

package info (click to toggle)
fastjet 3.0.6%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 9,404 kB
  • ctags: 3,770
  • sloc: cpp: 21,498; sh: 10,546; fortran: 673; makefile: 527; ansic: 131
file content (145 lines) | stat: -rw-r--r-- 5,207 bytes parent folder | download | duplicates (3)
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
/// \file
/// \page Examples FastJet examples
///
/// The FastJet examples have been organised by order of complexity,
/// starting by the simplest case and introducing features one after
/// another.
///   - \subpage Example01
///   - \subpage Example02
///   - \subpage Example03
///   - \subpage Example04
///   - \subpage Example05
///   - \subpage Example06
///   - \subpage Example07 (\subpage Example07old "old version")
///   - \subpage Example08
///   - \subpage Example09
///   - \subpage Example10
///   - \subpage Example11
///   - \subpage Example12 (\subpage Example12old "old version")
///   - \subpage Example13
///   - \subpage Example14

//STARTHEADER
// $Id: fastjet_example.cc 2684 2011-11-14 07:41:44Z soyez $
//
// Copyright (c) 2005-2011, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
//
//----------------------------------------------------------------------
// This file is part of FastJet.
//
//  FastJet is free software; you can redistribute it and/or modify
//  it under the terms of the GNU General Public License as published by
//  the Free Software Foundation; either version 2 of the License, or
//  (at your option) any later version.
//
//  The algorithms that underlie FastJet have required considerable
//  development and are described in hep-ph/0512210. If you use
//  FastJet as part of work towards a scientific publication, please
//  include a citation to the FastJet paper.
//
//  FastJet is distributed in the hope that it will be useful,
//  but WITHOUT ANY WARRANTY; without even the implied warranty of
//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
//  GNU General Public License for more details.
//
//  You should have received a copy of the GNU General Public License
//  along with FastJet. If not, see <http://www.gnu.org/licenses/>.
//----------------------------------------------------------------------
//ENDHEADER


//----------------------------------------------------------------------
// fastjet example program. 
//
// Compile it with: make fastjet_example
// run it with    : ./fastjet_example < data/single-event.dat
//
// People who are familiar with the ktjet package are encouraged to
// compare this file to the ktjet_example.cc program which does the
// same thing in the ktjet framework.
//----------------------------------------------------------------------
#include "fastjet/PseudoJet.hh"
#include "fastjet/ClusterSequence.hh"
#include<iostream> // needed for io
#include<sstream>  // needed for internal io
#include<vector> 
#include <cstdio>

using namespace std;

// a declaration of a function that pretty prints a list of jets
void print_jets (const vector<fastjet::PseudoJet> &);

/// an example program showing how to use fastjet
int main () {
  
  vector<fastjet::PseudoJet> input_particles;
  
  // Read in input particles
  double px, py , pz, E;
  while (cin >> px >> py >> pz >> E) {
    // create a fastjet::PseudoJet with these components and put it onto
    // back of the input_particles vector
    input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 
  }
  
  // create an object that represents your choice of jet algorithm and 
  // the associated parameters
  double Rparam = 1.0;
  fastjet::Strategy strategy = fastjet::Best;
  fastjet::RecombinationScheme recomb_scheme = fastjet::E_scheme;
  fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy);

  // run the jet clustering with the above jet definition
  fastjet::ClusterSequence clust_seq(input_particles, jet_def);

  // tell the user what was done
  cout << "Ran " << jet_def.description() << endl;
  cout << "Strategy adopted by FastJet was "<<
       clust_seq.strategy_string()<<endl<<endl;

  // extract the inclusive jets with pt > 5 GeV
  double ptmin = 5.0;
  vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);

  // print them out
  cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
  cout << "---------------------------------------\n";
  print_jets(inclusive_jets);
  cout << endl;

  // extract the exclusive jets with dcut = 25 GeV^2 
  double dcut = 25.0;
  vector<fastjet::PseudoJet> exclusive_jets = clust_seq.exclusive_jets(dcut);

  // print them out
  cout << "Printing exclusive jets with dcut = "<< dcut<<" GeV^2\n";
  cout << "--------------------------------------------\n";
  print_jets(exclusive_jets);


}


//----------------------------------------------------------------------
/// a function that pretty prints a list of jets
void print_jets (const vector<fastjet::PseudoJet> & jets) {

  // sort jets into increasing pt
  vector<fastjet::PseudoJet> sorted_jets = sorted_by_pt(jets);  

  // label the columns
  printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", 
	 "phi", "pt", "n constituents");
  
  // print out the details for each jet
  for (unsigned int i = 0; i < sorted_jets.size(); i++) {
    // the following is not super efficient since it creates an
    // intermediate constituents vector
    int n_constituents = sorted_jets[i].constituents().size();
    printf("%5u %15.8f %15.8f %15.8f %8u\n",
	   i, sorted_jets[i].rap(), sorted_jets[i].phi(),
	   sorted_jets[i].perp(), n_constituents);
  }

}