-
Notifications
You must be signed in to change notification settings - Fork 23
/
Copy pathrunSimpleJetAnalysis.cc
161 lines (117 loc) · 5.7 KB
/
runSimpleJetAnalysis.cc
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
#include <iostream>
#include <chrono>
#include "TFile.h"
#include "TTree.h"
#include "fastjet/PseudoJet.hh"
#include "fastjet/ClusterSequenceArea.hh"
#include "include/ProgressBar.h"
#include "PU14/EventMixer.hh"
#include "PU14/CmdLine.hh"
#include "PU14/PU14.hh"
#include "include/extraInfo.hh"
#include "include/jetCollection.hh"
#include "include/softDropGroomer.hh"
#include "include/treeWriter.hh"
#include "include/jetMatcher.hh"
#include "include/Angularity.hh"
using namespace std;
using namespace fastjet;
// ./runSimpleJetAnalysis -hard /Users/mverweij/mnt/eos/project/j/jetquenching/JetWorkshop2017/samples/pythia8/dijet120/PythiaEventsTune14PtHat120_0.pu14 -nev 10
int main (int argc, char ** argv) {
auto start_time = chrono::steady_clock::now();
CmdLine cmdline(argc,argv);
// inputs read from command line
int nEvent = cmdline.value<int>("-nev",1); // first argument: command line option; second argument: default value
//bool verbose = cmdline.present("-verbose");
cout << "will run on " << nEvent << " events" << endl;
// Uncomment to silence fastjet banner
ClusterSequence::set_fastjet_banner_stream(NULL);
//to write info to root tree
treeWriter trw("jetTree");
//Jet definition
double R = 0.4;
double ghostRapMax = 6.0;
double ghost_area = 0.005;
int active_area_repeats = 1;
GhostedAreaSpec ghost_spec(ghostRapMax, active_area_repeats, ghost_area);
AreaDefinition area_def = AreaDefinition(active_area,ghost_spec);
JetDefinition jet_def(antikt_algorithm, R);
double jetRapMax = 3.0;
Selector jet_selector = SelectorAbsRapMax(jetRapMax);
Angularity width(1.,1.,R);
Angularity pTD(0.,2.,R);
ProgressBar Bar(cout, nEvent);
Bar.SetStyle(-1);
EventMixer mixer(&cmdline); //the mixing machinery from PU14 workshop
// loop over events
int iev = 0;
unsigned int entryDiv = (nEvent > 200) ? nEvent / 200 : 1;
while ( mixer.next_event() && iev < nEvent )
{
// increment event number
iev++;
Bar.Update(iev);
Bar.PrintWithMod(entryDiv);
vector<PseudoJet> particlesMergedAll = mixer.particles();
vector<double> eventWeight;
eventWeight.push_back(mixer.hard_weight());
eventWeight.push_back(mixer.pu_weight());
// extract hard partons that initiated the jets
fastjet::Selector parton_selector = SelectorVertexNumber(-1);
vector<PseudoJet> partons = parton_selector(particlesMergedAll);
// select final state particles from hard event only
fastjet::Selector sig_selector = SelectorVertexNumber(0);
vector<PseudoJet> particlesSig = sig_selector(particlesMergedAll);
// select final state particles from background event only
fastjet::Selector bkg_selector = SelectorVertexNumber(1);
vector<PseudoJet> particlesBkg = bkg_selector(particlesMergedAll);
vector<PseudoJet> particlesMerged = particlesBkg;
particlesMerged.insert( particlesMerged.end(), particlesSig.begin(), particlesSig.end() );
//vector<PseudoJet> particlesBkg, particlesSig;
//SelectorIsHard().sift(particlesMerged, particlesSig, particlesBkg); // this sifts the full event into two vectors of PseudoJet, one for the hard event, one for the underlying event
//---------------------------------------------------------------------------
// jet clustering
//---------------------------------------------------------------------------
fastjet::ClusterSequenceArea csSig(particlesSig, jet_def, area_def);
jetCollection jetCollectionSig(sorted_by_pt(jet_selector(csSig.inclusive_jets(10.))));
//calculate some angularities
vector<double> widthSig; widthSig.reserve(jetCollectionSig.getJet().size());
vector<double> pTDSig; pTDSig.reserve(jetCollectionSig.getJet().size());
for(PseudoJet jet : jetCollectionSig.getJet()) {
widthSig.push_back(width.result(jet));
pTDSig.push_back(pTD.result(jet));
}
jetCollectionSig.addVector("widthSig", widthSig);
jetCollectionSig.addVector("pTDSig", pTDSig);
//---------------------------------------------------------------------------
// Groom the jets
//---------------------------------------------------------------------------
//SoftDrop grooming classic for signal jets (zcut=0.1, beta=0)
softDropGroomer sdgSigBeta00Z01(0.1, 0.0, R);
jetCollection jetCollectionSigSDBeta00Z01(sdgSigBeta00Z01.doGrooming(jetCollectionSig));
jetCollectionSigSDBeta00Z01.addVector("zgSigSDBeta00Z01", sdgSigBeta00Z01.getZgs());
jetCollectionSigSDBeta00Z01.addVector("ndropSigSDBeta00Z01", sdgSigBeta00Z01.getNDroppedSubjets());
jetCollectionSigSDBeta00Z01.addVector("dr12SigSDBeta00Z01", sdgSigBeta00Z01.getDR12());
//---------------------------------------------------------------------------
// write tree
//---------------------------------------------------------------------------
//Give variable we want to write out to treeWriter.
//Only vectors of the types 'jetCollection', and 'double', 'int', 'PseudoJet' are supported
trw.addCollection("eventWeight", eventWeight);
trw.addPartonCollection("partons", partons);
trw.addCollection("sigJet", jetCollectionSig);
trw.addCollection("sigJetSDBeta00Z01", jetCollectionSigSDBeta00Z01);
trw.fillTree();
}//event loop
Bar.Update(nEvent);
Bar.Print();
Bar.PrintLine();
TTree *trOut = trw.getTree();
TFile *fout = new TFile(cmdline.value<string>("-output", "JetToyHIResultSimpleJetAnalysis.root").c_str(), "RECREATE");
trOut->Write();
fout->Write();
fout->Close();
double time_in_seconds = chrono::duration_cast<chrono::milliseconds>
(chrono::steady_clock::now() - start_time).count() / 1000.0;
cout << "runFromFile: " << time_in_seconds << endl;
}