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