12 using namespace Pythia8;
15 #include "fastjet/PseudoJet.hh"
16 #include "fastjet/ClusterSequence.hh"
17 #include "fastjet/CDFMidPointPlugin.hh"
18 #include "fastjet/CDFJetCluPlugin.hh"
19 #include "fastjet/D0RunIIConePlugin.hh"
26 double pTfirstJet(
const Event& event,
int nJetMin,
double Rparam) {
28 double yPartonMax = 4.;
31 fastjet::Strategy strategy = fastjet::Best;
32 fastjet::RecombinationScheme recombScheme = fastjet::E_scheme;
33 fastjet::JetDefinition *jetDef = NULL;
35 if(event[3].colType() != 0 || event[4].colType() != 0)
36 jetDef =
new fastjet::JetDefinition(fastjet::kt_algorithm, Rparam,
37 recombScheme, strategy);
40 jetDef =
new fastjet::JetDefinition(fastjet::ee_kt_algorithm,
41 recombScheme, strategy);
43 std::vector <fastjet::PseudoJet> fjInputs;
48 for (
int i = 0; i <
event.size(); ++i) {
50 if ( !event[i].isFinal()
51 ||
event[i].isLepton()
52 ||
event[i].id() == 23
53 || abs(event[i].
id()) == 24
54 || abs(event[i].y()) > yPartonMax)
58 fjInputs.push_back( fastjet::PseudoJet (event[i].px(),
59 event[i].py(), event[i].pz(),event[i].e() ) );
63 if (
int(fjInputs.size()) == 0) {
69 fastjet::ClusterSequence clustSeq(fjInputs, *jetDef);
71 double pTFirst = sqrt(clustSeq.exclusive_dmerge_max(nJetMin-1));
95 virtual double tmsDefinition(
const Event& event);
98 double myKTdurham(
const Particle& RadAfterBranch,
99 const Particle& EmtAfterBranch,
int Type,
double D );
106 MyMergingHooks::MyMergingHooks() {}
109 MyMergingHooks::~MyMergingHooks() {}
115 double MyMergingHooks::tmsDefinition(
const Event& event){
119 int nFinalColoured = 0;
121 for(
int i=0; i <
event.size(); ++i) {
122 if(event[i].isFinal()){
123 if(event[i].
id() != 23 && abs(event[i].
id()) != 24)
125 if( event[i].colType() != 0)
131 int nLeptons = nHardOutLeptons();
132 int nQuarks = nHardOutPartons();
133 int nResNow = nResInCurrent();
136 if(nFinalNow - ( (nLeptons+nQuarks)/2 - nResNow)*2 != nFinalColoured){
140 if(nFinalNow != nFinalColoured)
return 0.;
144 int nMPI = infoPtr->nMPI();
145 if(nMPI > 1)
return 0.;
151 vector <int> FinalPartPos;
152 FinalPartPos.clear();
154 for (
int i=0; i <
event.size(); ++i)
155 if(event[i].isFinal() &&
event[i].colType() != 0)
156 FinalPartPos.push_back(i);
160 int type = (
event[3].colType() == 0 &&
event[4].colType() == 0) ? 1 : kTtype;
162 double ktmin =
event[0].e();
163 for(
int i=0; i < int(FinalPartPos.size()); ++i){
166 if(type == -1 || type == -2) {
167 double temp =
event[FinalPartPos[i]].pT();
168 kt12 = min(kt12, temp);
171 for(
int j=i+1; j < int(FinalPartPos.size()); ++j) {
172 double temp = kTdurham( event[FinalPartPos[i]], event[FinalPartPos[j]],
174 kt12 = min(kt12, temp);
177 ktmin = min(ktmin,kt12);
189 double MyMergingHooks::myKTdurham(
const Particle& RadAfterBranch,
190 const Particle& EmtAfterBranch,
int Type,
double D ){
195 Vec4 jet1 = RadAfterBranch.p();
196 Vec4 jet2 = EmtAfterBranch.p();
202 if (jet1.pAbs()*jet2.pAbs() <=0.) costh = 1.;
204 costh = costheta(jet1,jet2);
207 ktdur = 2.0*min( pow(jet1.e(),2) , (pow(jet2.e(),2)) )*(1.0 - costh);
208 }
else if( Type == -1 ){
210 double eta1 = 0.5*log( (jet1.e() + jet1.pz()) / (jet1.e() - jet1.pz()) );
211 double eta2 = 0.5*log( (jet2.e() + jet2.pz()) / (jet2.e() - jet2.pz()) );
213 double pt1 = sqrt( pow(jet1.px(),2) + pow(jet1.py(),2) );
214 double pt2 = sqrt( pow(jet2.px(),2) + pow(jet2.py(),2) );
215 double cosdPhi = ( jet1.px()*jet2.px() + jet1.py()*jet2.py() ) / (pt1*pt2);
216 double dPhi = acos( cosdPhi );
218 ktdur = min( pow(pt1,2),pow(pt2,2) )
219 * ( pow(eta1-eta2,2) + pow(dPhi,2) ) / pow(D,2);
220 }
else if( Type == -2 ){
222 double eta1 = 0.5*log( (jet1.e() + jet1.pz()) / (jet1.e() - jet1.pz()) );
223 double eta2 = 0.5*log( (jet2.e() + jet2.pz()) / (jet2.e() - jet2.pz()) );
224 double coshdEta = cosh( eta1 - eta2 );
226 double pt1 = sqrt( pow(jet1.px(),2) + pow(jet1.py(),2) );
227 double pt2 = sqrt( pow(jet2.px(),2) + pow(jet2.py(),2) );
228 double cosdPhi = ( jet1.px()*jet2.px() + jet1.py()*jet2.py() ) / (pt1*pt2);
230 ktdur = 2.0*min( pow(pt1,2),pow(pt2,2) )
231 * ( coshdEta - cosdPhi ) / pow(D,2);
243 int main(
int argc,
char* argv[] ){
247 cerr <<
" Unexpected number of command-line arguments. \n You are"
248 <<
" expected to provide the arguments \n"
249 <<
" 1. Input file for settings \n"
250 <<
" 2. Full name of the input LHE file (with path) \n"
251 <<
" 3. Path for output histogram files \n"
252 <<
" Program stopped. " << endl;
262 pythia.readFile(argv[1]);
263 string iPath = string(argv[2]);
264 string oPath = string(argv[3]);
267 int nEvent = pythia.mode(
"Main:numberOfEvents");
271 pythia.setMergingHooksPtr( myMergingHooks );
274 pythia.settings.forceParm(
"SpaceShower:pT0Ref",0.);
277 Hist histPTFirst(
"pT of first jet",100,0.,100.);
278 Hist histPTSecond(
"pT of second jet",100,0.,100.);
281 pythia.init(iPath,
false);
283 if(pythia.flag(
"Main:showChangedSettings")) {
284 pythia.settings.listChanged();
288 for(
int iEvent=0; iEvent<nEvent; ++iEvent ){
291 if( ! pythia.next())
continue;
294 double weight = pythia.info.mergingWeight();
297 double pTfirst = pTfirstJet(pythia.event,1, 0.4);
298 double pTsecnd = pTfirstJet(pythia.event,2, 0.4);
299 histPTFirst.fill( pTfirst, weight);
300 histPTSecond.fill( pTsecnd, weight);
302 if(iEvent%1000 == 0) cout << iEvent << endl;
311 * pythia.info.sigmaGen()
312 * 1./ double(nEvent);
315 histPTSecond *= norm;
318 string jetsInLHEF = iPath.substr(iPath.size()-5, iPath.size());
319 jetsInLHEF = jetsInLHEF.substr(0, jetsInLHEF.size()-4);
328 suffix << jetsInLHEF <<
"_wv.dat";
331 write.open( (
char*)(oPath +
"PTjet1_" + suffix.str()).c_str());
332 histPTFirst.table(write);
335 write.open( (
char*)(oPath +
"PTjet2_" + suffix.str()).c_str());
336 histPTSecond.table(write);
340 delete myMergingHooks;