ATLAS Offline Software
Loading...
Searching...
No Matches
Prophecy4fMerger.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4#include "Prophecy4fMerger.h"
5#include "TLorentzRotation.h"
6#include "TLorentzVector.h"
7#include "TVector3.h"
8/* Utils */
9#include "Pythia8/LHEF3.h"
10namespace LHEF {
11using Reader = Pythia8::Reader;
12using Writer = Pythia8::Writer;
13}
14#include <cmath>
15#include <iostream>
16#include <algorithm>
18
19using namespace asg::msgUserCode;
20
23
26
27void Prophecy4fMerger::setRandomSeed(unsigned long long seed)
28{
29 m_rand.SetSeed(seed);
30 ANA_MSG_INFO("Prophecy4fMerger::setRandomSeed: set random seed to " << seed << ".");
31}
32
33int Prophecy4fMerger::alulb4(double *ps, double *pi, double *pf){
34
35 // transform pi rest to lab frame of ps - output pf
36 // For description, see LOREN4 doc in CERNLIB
37
38 if(ps[3] == ps[4]){
39 pf[0] = pi[0];
40 pf[1] = pi[1];
41 pf[2] = pi[2];
42 pf[3] = pi[3];
43 }
44 else{
45 const double pf4 = (pi[0] * ps[0] + pi[1] * ps[1] + pi[2] * ps[2] + pi[3] * ps[3]) / ps[4];
46 const double fn = (pf4 + pi[3]) / (ps[3] + ps[4]);
47 pf[0] = pi[0] + fn * ps[0];
48 pf[1] = pi[1] + fn * ps[1];
49 pf[2] = pi[2] + fn * ps[2];
50 pf[3] = pf4;
51 }
52 return 0;
53
54}
55
56int Prophecy4fMerger::alulf4(double *ps, double *pi, double *pf){
57
58 // transform pi lab to rest frame of ps - output pf
59 // For description, see LOREN4 doc in CERNLIB
60
61 if(ps[3] == ps[4]){
62 pf[0] = pi[0];
63 pf[1] = pi[1];
64 pf[2] = pi[2];
65 pf[3] = pi[3];
66 }
67 else{
68 const double pf4 = (pi[3] * ps[3] - pi[2] * ps[2] - pi[1] * ps[1] - pi[0] * ps[0]) / ps[4];
69 const double fn = (pf4 + pi[3]) / (ps[3] + ps[4]);
70 pf[0] = pi[0] - fn * ps[0];
71 pf[1] = pi[1] - fn * ps[1];
72 pf[2] = pi[2] - fn * ps[2];
73 pf[3] = pf4;
74 }
75 return 0;
76
77}
78
79int Prophecy4fMerger::alulob(double *ps, double *pi, double *pf){
80
81 //extern int alulb4(double *, double *, double *);
82
83 alulb4(ps, pi, pf);
84 pf[4] = pi[4];
85 return 0;
86
87}
88
89int Prophecy4fMerger::alulof(double *ps, double *pi, double *pf){
90
91 //extern int alulf4(double *, double *, double *);
92
93 alulf4(ps, pi, pf);
94 pf[4] = pi[4];
95
96 return 0;
97
98}
99
100double Prophecy4fMerger::alupcm(double em0, double em1, double em2){
101
102 // calculate the momentum of particles 1 and 2 in the rest frame of 0 given the three masses
103
104 double ret_val;
105
106 //ems = abs(em1 + em2);
107 //emd = abs(em1 - em2);
108 double ems = (em1 + em2);
109 if(ems < 0){
110 ems=-(em1 + em2);
111 }
112 double emd = (em1 - em2);
113 if(emd < 0){
114 emd = -(em1 - em2);
115 }
116
117 if(em0 < ems || em0 < emd){
118 ret_val = -1.f;
119 }
120 else if(em0 == ems || em0 == emd){
121 ret_val = 0.f;
122 }
123 else{
124 ret_val = sqrt((em0 + emd) * (em0 - emd) * (em0 + ems) * (em0 - ems)) * .5f / em0;
125 }
126 return ret_val;
127
128}
129
130int Prophecy4fMerger::rescms(double *p, double *p1, double *p2, double m1, double m2){
131
132 // p1 and p2 are the four-vector decay products of p. Here we change the masses of p1 and p2 to be
133 // m1 and m2. Note that for Prophecy, the masses in p1 and p2 are zero.
134
135 double po1[5], po2[5], pcm, pcmo;
136
137 const double m = p[4];
138 const double mo1 = p1[4];
139 const double mo2 = p2[4];
140
141 // transform p1 to po1 to be in rest frame of p, sampe for p2
142 alulof(p, p1, po1);
143 alulof(p, p2, po2);
144
145 // Calculate the momentum of p1/p2 in the rest frame of p (po1/po2) for the original mass and new
146 // masses, and then rescale the rest frame momenta po1 and po2. (Note: eo1 = eo2 in cms of p, so
147 // po1 and po2 in cms of p change according to e^2 - m^2 as m changes.)
148 if (std::max(mo1,mo2) < 1e-6){
149 pcmo = m / 2.;
150 }
151 else{
152 pcmo = alupcm(m, mo1, mo2);
153 }
154 pcm = alupcm(m, m1, m2);
155
156 //rescale the cms momenta, po1 and po2, to account for the new masses used
157 for (int il = 0; il < 4; il++) {
158 po1[il] = pcm / pcmo * po1[il];
159 po2[il] = pcm / pcmo * po2[il];
160 }
161
162 //set energy and mass
163 po1[3] = sqrt(pcm*pcm + m1*m1);
164 po2[3] = sqrt(pcm*pcm + m2*m2);
165 po1[4] = m1;
166 po2[4] = m2;
167
168 //boost back
169 alulob(p, po1, p1);
170 alulob(p, po2, p2);
171
172 return 0;
173
174}
175
176void Prophecy4fMerger::setIO(const std::string& powheg,
177 const std::string& prophecy4e,
178 const std::string& prophecy4mu,
179 const std::string& prophecy2e2mu,
180 const std::string& outlhe,
181 bool debug)
182{
183 m_inPowheg = powheg;
184 m_inProphecy4e = prophecy4e;
185 m_inProphecy4mu = prophecy4mu;
186 m_inProphecy2e2mu = prophecy2e2mu;
187 m_outLHE = outlhe;
188 m_phEvent = false;
189 m_debug = debug;
190
191}
192
194
195 print(" Opening Powheg LHE file ... " + m_inPowheg);
196 if( !fileExists(m_inPowheg) ){
197 ANA_MSG_ERROR("Input Powheg LHE not Found!"
198 << " Aborting ... ");
199 return;
200 }
201 LHEF::Reader readH( m_inPowheg.c_str() );
202
203 print(" Opening Prophecy4f LHE file for 4e ... " + m_inProphecy4e);
205 ANA_MSG_INFO( "Input Prophecy4f for 4e LHE not Found! " << m_inProphecy4e
206 << " Aborting ... " );
207 return;
208 }
209 LHEF::Reader read4e( m_inProphecy4e.c_str() );
210
211 print(" Opening Prophecy4f LHE file for 4mu ... " + m_inProphecy4mu);
213 ANA_MSG_ERROR("Input Prophecy4f for 4mu LHE not Found! " << m_inProphecy4mu
214 << " Aborting ... ");
215 return;
216 }
217 LHEF::Reader read4mu( m_inProphecy4mu.c_str() );
218
219 print(" Opening Prophecy4f for 2e2mu LHE file ... " + m_inProphecy2e2mu);
221 ANA_MSG_ERROR("Input Prophecy4f for 2e2mu LHE not Found! " << m_inProphecy2e2mu
222 << " Aborting ... ");
223 return;
224 }
225 LHEF::Reader read2e2mu( m_inProphecy2e2mu.c_str() );
226
227 print(" Opening Out LHE file ... " + m_outLHE);
228 LHEF::Writer writeLHE( m_outLHE.c_str() );
229
230 writeLHE.headerBlock() << readH.headerBlock;
231 writeLHE.initComments() << readH.initComments;
232 writeLHE.heprup = readH.heprup;
233 writeLHE.init();
234
235 long neve = 0;
236 LHEF::Reader* read4f = 0;
237
238
239 // Read in each Higgs event and one Prophecy4f event, merge them and write out the merged event
240 // randomly chose between 4mu, 4e, 2mu2e and 2e2mu. Note that 2mu2e is always taken as 2e2mu,
241 // since it is generated this way.
242
243 while( readH.readEvent() ){
244 // randomly choose between 4mu, 4e, 2mu2e and 2e2mu
245 int ID_l[4] = {-999, -999, -999, -999};
246 int decayModeType = 4.0 * m_rand.Rndm();
247
248 if (decayModeType == 0){
249 //4e
250 read4f = &read4e;
251 for(int i = 0; i < 4; i++) ID_l[i] = m_electronID;
252 }
253 else if (decayModeType == 1 || decayModeType == 2){
254 //2e2mu
255 read4f = &read2e2mu;
256 for(int i = 0; i < 2; i++) ID_l[i] = m_electronID;
257 for(int i = 2; i < 4; i++) ID_l[i] = m_muonID;
258 }
259 else {
260 //4mu
261 read4f = &read4mu;
262 for(int i = 0; i < 4; i++) ID_l[i] = m_muonID;
263 }
264
265 if( !read4f->readEvent() ){
266 ANA_MSG_INFO("Still Powheg events but no more Prophecy4f ");
267 ANA_MSG_INFO("events! Writing out LHE file ... Events processed so far " << neve );
268 break;
269 }
270
271 ++neve;
272 if(neve%1000==0){
273 print("Events processed so far ", neve);
274 }
275 if( readH.outsideBlock.length() )
276 ANA_MSG_INFO( readH.outsideBlock );
277
278 int nup_org=readH.hepeup.NUP;
279
280 for(int nup=0; nup<nup_org; nup++){
281 if(readH.hepeup.IDUP[nup]==m_higgsID)
282 readH.hepeup.ISTUP[nup]=2; // set Higgs status code to 2
283 }
284 writeLHE.hepeup = readH.hepeup;
285
286 double P4_Z[5][2]{}; // momentum and mass of Zi i=1,2 set to 0 P4_Z[5] = (px, py, pz, E, mZ)
287 double P4_l[5][5]{};
288
289 double Pph[5] = {-999., -999., -999., -999., -999.};
290 double daughter2[5][5]{};
291
292 for(int i=0; i<nup_org; i++){
293 if(readH.hepeup.IDUP[i]==m_higgsID) {
294 TLorentzVector higgs( readH.hepeup.PUP[i][0],
295 readH.hepeup.PUP[i][1],
296 readH.hepeup.PUP[i][2],
297 readH.hepeup.PUP[i][3] );
298
299 TVector3 boostvec = higgs.BoostVector();
300 TLorentzRotation boost(-boostvec);
301 TLorentzRotation backboost(boostvec);
302
303 /* boost to Higgs rest frame */
304 TLorentzVector r_higgs(higgs);
305 r_higgs*=boost;
306
307 /* Here add daughter particles to Higgs in rest
308 boost back with backboost */
309 writeLHE.hepeup.resize(nup_org+read4f->hepeup.NUP);
310
311 TLorentzVector sum_daughter;
312 TLorentzVector sum_daughter_rest;
313 TLorentzVector sum_daughter_rest_init;
314
315 for(int k=0; k<read4f->hepeup.NUP; k++){
316 TLorentzVector daughter(read4f->hepeup.PUP[k][0],
317 read4f->hepeup.PUP[k][1],
318 read4f->hepeup.PUP[k][2],
319 read4f->hepeup.PUP[k][3]);
320 sum_daughter_rest_init+=daughter;
321 }
322
323 for(int j=0; j<read4f->hepeup.NUP; j++){
324 TLorentzVector daughter(read4f->hepeup.PUP[j][0],
325 read4f->hepeup.PUP[j][1],
326 read4f->hepeup.PUP[j][2],
327 read4f->hepeup.PUP[j][3]);
328 sum_daughter_rest+=daughter;
329
330 daughter *= backboost;
331 daughter *= higgs.M()/sum_daughter_rest_init.M();
332
333 /* POWHEG event with Higgs off-mass shell */
334 m_phEvent=isPHevent(higgs,sum_daughter_rest_init);
335 if( m_phEvent ){
336 break;
337 }
338
339 sum_daughter+=daughter;
340
341 if( std::abs(read4f->hepeup.IDUP[j])==m_electronID ||
342 std::abs(read4f->hepeup.IDUP[j])==m_muonID ){
343
344 //ID_l[j] = read4f->hepeup.IDUP[j];
345 P4_l[0][j] = daughter.Px();
346 P4_l[1][j] = daughter.Py();
347 P4_l[2][j] = daughter.Pz();
348 P4_l[3][j] = daughter.E();
349
350 }
351 else if( std::abs(read4f->hepeup.IDUP[j])==m_photonID ){
352 Pph[0] = daughter.Px();
353 Pph[1] = daughter.Py();
354 Pph[2] = daughter.Pz();
355 Pph[3] = daughter.E();
356 Pph[4] =0.0;
357 }
358 }
359
360 for(int n=0; n<4; n++){
361 if( std::abs(ID_l[n])==m_electronID ||
362 std::abs(ID_l[n])==m_muonID ) {
363
364 // Z1 -> l1 + l2 and Z2 -> l3 + l4
365 P4_Z[n][0] = P4_l[n][0] + P4_l[n][1];
366 P4_Z[n][1] = P4_l[n][2] + P4_l[n][3];
367 }
368 }
369
370 double mass1, mass2;
371 P4_Z[4][0] = sqrt(P4_Z[3][0] * P4_Z[3][0] - P4_Z[0][0] * P4_Z[0][0] - P4_Z[1][0] * P4_Z[1][0] - P4_Z[2][0] * P4_Z[2][0]);
372 P4_Z[4][1] = sqrt(P4_Z[3][1] * P4_Z[3][1] - P4_Z[0][1] * P4_Z[0][1] - P4_Z[1][1] * P4_Z[1][1] - P4_Z[2][1] * P4_Z[2][1]);
373 mass1 = setParticleMass(ID_l[0]);
374 mass2 = setParticleMass(ID_l[2]);
375
376 double p1[5], p2[5], pZ1[5], p3[5], p4[5], pZ2[5];
377 for(int u=0; u<5; u++){
378
379 pZ1[u] = P4_Z[u][0];
380 p1[u] = P4_l[u][0];
381 p2[u] = P4_l[u][1];
382 pZ2[u] = P4_Z[u][1];
383 p3[u] = P4_l[u][2];
384 p4[u] = P4_l[u][3];
385
386 }
387
388 p1[4] = 0.0;
389 p2[4] = 0.0;
390 p3[4] = 0.0;
391 p4[4] = 0.0;
392
393 // Give masses to the leptons using Alpgen/FHerwig routines
394 // Note: prophecy leptons are generated massless
395 rescms(pZ1, p1, p2, mass1, mass1);
396 rescms(pZ2, p3, p4, mass2, mass2);
397
398 for(int d =0; d<5; d++){
399 daughter2[0][d] = p1[d];
400 daughter2[1][d] = p2[d];
401 daughter2[2][d] = p3[d];
402 daughter2[3][d] = p4[d];
403 daughter2[4][d] = Pph[d];
404 }
405
406 // print out children
407 bool doDebug = false;
408 doDebug = true;
409 if (doDebug) {
410
411 ANA_MSG_DEBUG("Event: " << neve);
412 ANA_MSG_DEBUG(" Higgs and children, m, pt, eta, phi" );
413 ANA_MSG_DEBUG("Higgs " << higgs.M() << ", " << higgs.Pt() << "," << higgs.Eta() << "," << higgs.Phi());
414
415
416 TLorentzVector higgsFromChildren;
417 for(int d =0; d < read4f->hepeup.NUP; d++){
418 TLorentzVector child;
419 child.SetXYZM( daughter2[d][0],
420 daughter2[d][1],
421 daughter2[d][2],
422 daughter2[d][4] );
423
424 ANA_MSG_DEBUG("child " << d << " " << child.M() << ", " << child.Pt() << ", " << child.Eta() << ", " << child.Phi());
425 higgsFromChildren += child;
426 }
427 ANA_MSG_DEBUG("Higgs mass " << higgs.M() << ", from children " << higgsFromChildren.M() << ", diff " << higgs.M() - higgsFromChildren.M());
428 }
429 break;
430 }
431 }
432
433 /* Write LHE */
434 if(m_phEvent){
435 writeLHE.hepeup = readH.hepeup;
436 }
437 else{
438
439 for(int jp=0; jp<read4f->hepeup.NUP; jp++){
440 if (jp < 4){
441 writeLHE.hepeup.IDUP[nup_org+jp] = (jp % 2 == 0)? ID_l[jp] : -ID_l[jp];
442
443 }
444 else{
445 writeLHE.hepeup.IDUP[nup_org+jp] = read4f->hepeup.IDUP[jp];
446 }
447
448 writeLHE.hepeup.ISTUP[nup_org+jp]=read4f->hepeup.ISTUP[jp];
449 if( read4f->hepeup.ICOLUP[jp].first==0 &&
450 read4f->hepeup.ICOLUP[jp].second==0 ) {
451 writeLHE.hepeup.ICOLUP[nup_org+jp].first=0;
452 writeLHE.hepeup.ICOLUP[nup_org+jp].second=0;
453 }
454 else {
455 writeLHE.hepeup.ICOLUP[nup_org+jp].first=read4f->hepeup.ICOLUP[jp].first+100;
456 writeLHE.hepeup.ICOLUP[nup_org+jp].second=read4f->hepeup.ICOLUP[jp].second+100;
457 }
458 writeLHE.hepeup.MOTHUP[nup_org+jp].first=3;//i-1;
459 writeLHE.hepeup.MOTHUP[nup_org+jp].second=3;//i-1;
460 writeLHE.hepeup.PUP[nup_org+jp][0]=daughter2[jp][0];
461 writeLHE.hepeup.PUP[nup_org+jp][1]=daughter2[jp][1];
462 writeLHE.hepeup.PUP[nup_org+jp][2]=daughter2[jp][2];
463 writeLHE.hepeup.PUP[nup_org+jp][3]=daughter2[jp][3];
464 writeLHE.hepeup.PUP[nup_org+jp][4]=daughter2[jp][4];
465 writeLHE.hepeup.VTIMUP[nup_org+jp]=read4f->hepeup.VTIMUP[jp];
466 writeLHE.hepeup.SPINUP[nup_org+jp]=read4f->hepeup.SPINUP[jp];
467 }
468 }
469
470 // set powheg weights to negative, if prophecy weight is -1
471
472 // Get Prophecy4f weight
473 writeLHE.hepeup.XWGTUP *= read4f->hepeup.XWGTUP;
474 double prophecyWeight = read4f->hepeup.XWGTUP;
475 if ( writeLHE.hepeup.weights_compressed.size() > 0 ) {
476 for ( int i = 1, N = writeLHE.hepeup.weights_compressed.size(); i < N; ++i ) {
477 writeLHE.hepeup.weights_compressed[i] *= prophecyWeight;
478
479 }
480 }
481
482 writeLHE.eventComments() << readH.eventComments;
483 writeLHE.writeEvent();
484 }
485}
486
487bool Prophecy4fMerger::isPHevent(const TLorentzVector& higgs,
488 const TLorentzVector& sum_daugh_rest_init){
489
490 if( std::abs(higgs.M()-sum_daugh_rest_init.M())>m_deltaM ){
491 ANA_MSG_INFO("POWHEG event with Higgs off-mass shell: " << higgs.M() << " GeV");
492 return true;
493 }
494 return false;
495
496}
497
498
500
501 double mass = 0.0;
502
503 if( std::abs(id) == m_neutrinoEl ||
504 std::abs(id) == m_neutrinoMu ||
505 std::abs(id) == m_neutrinoTau ){
506 mass = 0.0;
507 }
508 else if( std::abs(id) == m_electronID ){
509 mass = m_electronMass;
510 }
511 else if( std::abs(id) == m_muonID ){
512 mass = m_muonMass;
513 }
514 else if( std::abs(id) == m_tauID ){
515 mass = m_tauMass;
516 }
517
518 return mass;
519
520}
521
522bool Prophecy4fMerger::fileExists(const std::string& filename){
523
524 std::ifstream ifile(filename.c_str());
525 if( ifile.good() ){
526 ifile.close();
527 return true;
528 }
529 ifile.close();
530 return false;
531
532}
533
534void Prophecy4fMerger::print(const std::string& field){
535 if(m_debug){
536 ANA_MSG_DEBUG(field);
537 }
538}
539
540void Prophecy4fMerger::print(const std::string& field,
541 int value){
542
543 if(m_debug){
544 ANA_MSG_DEBUG("DEBUG:: " << value << " " << field);
545 }
546
547}
548
macros for messaging and checking status codes
#define ANA_MSG_INFO(xmsg)
Macro printing info messages.
#define ANA_MSG_ERROR(xmsg)
Macro printing error messages.
#define ANA_MSG_DEBUG(xmsg)
Macro printing debug messages.
const bool debug
#define pi
double alupcm(double em0, double em1, double em2)
bool fileExists(const std::string &filename)
static const long m_higgsID
static constexpr double m_deltaM
static constexpr double m_muonMass
int alulof(double *ps, double *pi, double *pf)
static const long m_neutrinoTau
std::string m_inProphecy4mu
std::string m_inProphecy2e2mu
static const long m_neutrinoEl
int alulob(double *ps, double *pi, double *pf)
void setIO(const std::string &powheg, const std::string &prophecy4e, const std::string &prophecy4mu, const std::string &prophecy2e2mu, const std::string &outlhe, bool debug)
std::string m_inProphecy4e
static const long m_neutrinoMu
int alulf4(double *ps, double *pi, double *pf)
static constexpr double m_electronMass
int alulb4(double *ps, double *pi, double *pf)
static constexpr double m_tauMass
static const long m_muonID
bool isPHevent(const TLorentzVector &higgs, const TLorentzVector &sum_daugh_rest_init)
static const long m_tauID
void setRandomSeed(unsigned long long seed)
int rescms(double *p, double *p1, double *p2, double m1, double m2)
void print(const std::string &field)
static const long m_photonID
static const long m_electronID
double setParticleMass(int id) const
Pythia8::Writer Writer
Pythia8::Reader Reader