ATLAS Offline Software
Loading...
Searching...
No Matches
EvtPythiaEngine.cxx
Go to the documentation of this file.
1
2/***********************************************************************
3* Copyright 1998-2023 CERN for the benefit of the EvtGen authors *
4* *
5* This file is part of EvtGen. *
6* *
7* EvtGen is free software: you can redistribute it and/or modify *
8* it under the terms of the GNU General Public License as published by *
9* the Free Software Foundation, either version 3 of the License, or *
10* (at your option) any later version. *
11* *
12* EvtGen is distributed in the hope that it will be useful, *
13* but WITHOUT ANY WARRANTY; without even the implied warranty of *
14* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15* GNU General Public License for more details. *
16* *
17* You should have received a copy of the GNU General Public License *
18* along with EvtGen. If not, see <https://www.gnu.org/licenses/>. *
19***********************************************************************/
20
21
22#include "EvtGen_i/EvtGenExternal/EvtPythiaEngine.hh"
23
24#include "EvtGenBase/EvtDecayTable.hh"
25#include "EvtGenBase/EvtExtGeneratorCommandsTable.hh"
26#include "EvtGenBase/EvtPDL.hh"
27#include "EvtGenBase/EvtParticleFactory.hh"
28#include "EvtGenBase/EvtReport.hh"
29#include "EvtGenBase/EvtSpinType.hh"
30
31#include "EvtGen_i/EvtGenExternal/EvtPythia6CommandConverter.hh"
32
33#include "Pythia8/Event.h"
34#include "Pythia8/ParticleData.h"
35
36#include <cmath>
37#include <iostream>
38#include <sstream>
39
40#if PYTHIA_VERSION_INTEGER < 8304
41typedef Pythia8::ParticleDataEntry* ParticleDataEntryPtr;
42#else
43typedef Pythia8::ParticleDataEntryPtr ParticleDataEntryPtr;
44#endif
45
46using std::endl;
47
48EvtPythiaEngine::EvtPythiaEngine( std::string xmlDir, bool convertPhysCodes,
49 bool useEvtGenRandom )
50{
51 // Create two Pythia generators. One will be for generic
52 // Pythia decays in the decay.dec file. The other one will be to
53 // only decay aliased particles, which are in general "signal"
54 // decays different from those in the decay.dec file.
55 // Even though it is not possible to have two different particle
56 // versions in one Pythia generator, we can use two generators to
57 // get the required behaviour.
58
59 EvtGenReport( EVTGEN_INFO, "EvtGen" )
60 << "Creating generic Pythia generator" << endl;
61 m_genericPythiaGen = std::make_unique<Pythia8::Pythia>( xmlDir );
62
63 EvtGenReport( EVTGEN_INFO, "EvtGen" )
64 << "Creating alias Pythia generator" << endl;
65 m_aliasPythiaGen = std::make_unique<Pythia8::Pythia>( xmlDir, false );
66
67 m_thePythiaGenerator = 0;
68 m_daugPDGVector.clear();
69 m_daugP4Vector.clear();
70
71 m_convertPhysCodes = convertPhysCodes;
72
73 // Specify if we are going to use the random number generator (engine)
74 // from EvtGen for Pythia 8.
75 m_useEvtGenRandom = useEvtGenRandom;
76
77 m_evtgenRandom = std::make_shared<EvtPythiaRandom>();
78
79 m_initialised = false;
80}
81
82EvtPythiaEngine::~EvtPythiaEngine()
83{
84 m_thePythiaGenerator = nullptr;
85 this->clearDaughterVectors();
86 this->clearPythiaModeMap();
87}
88
89void EvtPythiaEngine::clearDaughterVectors()
90{
91 m_daugPDGVector.clear();
92 m_daugP4Vector.clear();
93}
94
95void EvtPythiaEngine::clearPythiaModeMap()
96{
97 PythiaModeMap::iterator iter;
98 for ( iter = m_pythiaModeMap.begin(); iter != m_pythiaModeMap.end(); ++iter ) {
99 std::vector<int> modeVector = iter->second;
100 modeVector.clear();
101 }
102
103 m_pythiaModeMap.clear();
104}
105
106void EvtPythiaEngine::initialise()
107{
108 if ( m_initialised ) {
109 return;
110 }
111
112 this->clearPythiaModeMap();
113
114 this->updateParticleLists();
115
116 // Hadron-level processes only (hadronized, string fragmentation and secondary decays).
117 // We do not want to generate the full pp or e+e- event structure etc..
118 m_genericPythiaGen->readString( "ProcessLevel:all = off" );
119 m_aliasPythiaGen->readString( "ProcessLevel:all = off" );
120
121 // Turn off Pythia warnings, e.g. changes to particle properties
122 m_genericPythiaGen->readString( "Print:quiet = on" );
123 m_aliasPythiaGen->readString( "Print:quiet = on" );
124
125 // Apply any other physics (or special particle) requirements/cuts etc..
126 this->updatePhysicsParameters();
127
128 // Set the random number generator
129 if ( m_useEvtGenRandom == true ) {
130#if PYTHIA_VERSION_INTEGER >= 8310
131 m_genericPythiaGen->setRndmEnginePtr( m_evtgenRandom );
132 m_aliasPythiaGen->setRndmEnginePtr( m_evtgenRandom );
133#else
134 m_genericPythiaGen->setRndmEnginePtr( m_evtgenRandom.get() );
135 m_aliasPythiaGen->setRndmEnginePtr( m_evtgenRandom.get() );
136#endif
137 }
138
139 m_genericPythiaGen->init();
140 m_aliasPythiaGen->init();
141
142 m_initialised = true;
143}
144
145bool EvtPythiaEngine::doDecay( EvtParticle* theParticle )
146{
147 // Store the mother particle within a Pythia8 Event object.
148 // Then do the hadron level decays.
149 // The EvtParticle must be a colour singlet (meson/baryon/lepton), i.e. not a gluon or quark
150
151 // We delete any daughters the particle may have, since we are asking Pythia
152 // to generate the decay anew. Also note that _any_ Pythia decay allowed for the particle
153 // will be generated and not the specific Pythia decay mode that EvtGen has already
154 // specified. This is necessary since we only want to initialise the Pythia decay table
155 // once; all Pythia branching fractions for a given mother particle are renormalised to sum to 1.0.
156 // In EvtGen decay.dec files, it may be the case that Pythia decays are only used
157 // for some of the particle decays (i.e. Pythia BF sum < 1.0). As we loop over many events,
158 // the total frequency for each Pythia decay mode will normalise correctly to what
159 // we wanted via the specifications made to the decay.dec file, even though event-by-event
160 // the EvtGen decay channel and the Pythia decay channel may be different.
161
162 if ( m_initialised == false ) {
163 this->initialise();
164 }
165
166 if ( theParticle == 0 ) {
167 EvtGenReport( EVTGEN_INFO, "EvtGen" )
168 << "Error in EvtPythiaEngine::doDecay. The mother particle is null. Not doing any Pythia decay."
169 << endl;
170 return false;
171 }
172
173 // Delete EvtParticle daughters if they already exist
174 if ( theParticle->getNDaug() != 0 ) {
175 bool keepChannel( false );
176 theParticle->deleteDaughters( keepChannel );
177 }
178
179 EvtId particleId = theParticle->getId();
180 int isAlias = particleId.isAlias();
181
182 // Choose the generator depending if we have an aliased (parent) particle or not
183 m_thePythiaGenerator = ( isAlias == 1 ? m_aliasPythiaGen.get()
184 : m_genericPythiaGen.get() );
185
186 // Need to use the reference to the Pythia8::Event object,
187 // otherwise it will just return a new empty, default event object.
188 Pythia8::Event& theEvent = m_thePythiaGenerator->event;
189 theEvent.reset();
190
191 // Initialise the event to be the particle rest frame
192 int PDGCode = EvtPDL::getStdHep( particleId );
193
194 int status( 1 );
195 int colour( 0 ), anticolour( 0 );
196 double px( 0.0 ), py( 0.0 ), pz( 0.0 );
197 double m0 = theParticle->mass();
198 double E = m0;
199
200 theEvent.append( PDGCode, status, colour, anticolour, px, py, pz, E, m0 );
201
202 // Generate the Pythia event
203 int iTrial( 0 );
204 bool generatedEvent( false );
205 for ( iTrial = 0; iTrial < 10; iTrial++ ) {
206 generatedEvent = m_thePythiaGenerator->next();
207 if ( generatedEvent ) {
208 break;
209 }
210 }
211
212 bool success( false );
213
214 if ( generatedEvent ) {
215 // Store the daughters for this particle from the Pythia decay tree.
216 // This is a recursive function that will continue looping through
217 // all available daughters until the first set of non-quark and non-gluon
218 // particles are encountered in the Pythia Event structure.
219
220 // First, clear up the internal vectors storing the daughter
221 // EvtId types and 4-momenta.
222 this->clearDaughterVectors();
223
224 // Now store the daughter info. Since this is a recursive function
225 // to loop through the full Pythia decay tree, we do not want to create
226 // EvtParticles here but in the next step.
227 this->storeDaughterInfo( theParticle, 1 );
228
229 // Now create the EvtParticle daughters of the (parent) particle.
230 // We need to use the EvtParticle::makeDaughters function
231 // owing to the way EvtParticle stores parent-daughter information.
232 this->createDaughterEvtParticles( theParticle );
233
234 //theParticle->printTree();
235 //theEvent.list(true, true);
236
237 success = true;
238 }
239
240 return success;
241}
242
243void EvtPythiaEngine::storeDaughterInfo( EvtParticle* theParticle, int startInt )
244{
245 Pythia8::Event& theEvent = m_thePythiaGenerator->event;
246
247 std::vector<int> daugList = theEvent.daughterList( startInt );
248
249 std::vector<int>::iterator daugIter;
250 for ( daugIter = daugList.begin(); daugIter != daugList.end(); ++daugIter ) {
251 int daugInt = *daugIter;
252
253 // Ask if the daughter is a quark or gluon. If so, recursively search again.
254 Pythia8::Particle& daugParticle = theEvent[daugInt];
255
256 if ( daugParticle.isQuark() || daugParticle.isGluon() ) {
257 // Recursively search for correct daughter type
258 this->storeDaughterInfo( theParticle, daugInt );
259
260 } else {
261 // We have a daughter that is not a quark nor gluon particle.
262 // Make sure we are not double counting particles, since several quarks
263 // and gluons make one particle.
264 // Set the status flag for any "new" particle to say that we have stored it.
265 // Use status flag = 1000 (within the user allowed range for Pythia codes).
266
267 // Check that the status flag for the particle is not equal to 1000
268 int statusCode = daugParticle.status();
269 if ( statusCode != 1000 ) {
270 int daugPDGInt = daugParticle.id();
271
272 double px = daugParticle.px();
273 double py = daugParticle.py();
274 double pz = daugParticle.pz();
275 double E = daugParticle.e();
276 EvtVector4R daughterP4( E, px, py, pz );
277
278 // Now store the EvtId and 4-momentum in the internal vectors
279 m_daugPDGVector.push_back( daugPDGInt );
280 m_daugP4Vector.push_back( daughterP4 );
281
282 // Set the status flag for the Pythia particle to let us know
283 // that we have already considered it to avoid double counting.
284 daugParticle.status( 1000 );
285
286 } // Status code != 1000
287 }
288 }
289}
290
291void EvtPythiaEngine::createDaughterEvtParticles( EvtParticle* theParent )
292{
293 if ( theParent == 0 ) {
294 EvtGenReport( EVTGEN_INFO, "EvtGen" )
295 << "Error in EvtPythiaEngine::createDaughterEvtParticles. The parent is null"
296 << endl;
297 return;
298 }
299
300 // Get the list of Pythia decay modes defined for this particle id alias.
301 // It would be easier to just use the decay channel number that Pythia chose to use
302 // for the particle decay, but this is not accessible from the Pythia interface at present.
303
304 int nDaughters = m_daugPDGVector.size();
305 std::vector<EvtId> daugAliasIdVect( 0 );
306
307 EvtId particleId = theParent->getId();
308 // Check to see if we have an anti-particle. If we do, charge conjugate the particle id to get the
309 // Pythia "alias" we can compare with the defined (particle) Pythia modes.
310 int PDGId = EvtPDL::getStdHep( particleId );
311 int aliasInt = particleId.getAlias();
312 int pythiaAliasInt( aliasInt );
313
314 if ( PDGId < 0 ) {
315 // We have an anti-particle.
316 EvtId conjPartId = EvtPDL::chargeConj( particleId );
317 pythiaAliasInt = conjPartId.getAlias();
318 }
319
320 std::vector<int> pythiaModes = m_pythiaModeMap[pythiaAliasInt];
321
322 // Loop over all available Pythia decay modes and find the channel that matches
323 // the daughter ids. Set each daughter id to also use the alias integer.
324 // This will then convert the Pythia generated channel to the EvtGen alias defined one.
325
326 std::vector<int>::iterator modeIter;
327 bool gotMode( false );
328
329 for ( modeIter = pythiaModes.begin(); modeIter != pythiaModes.end();
330 ++modeIter ) {
331 // Stop the loop if we have the right decay mode channel
332 if ( gotMode ) {
333 break;
334 }
335
336 int pythiaModeInt = *modeIter;
337
338 EvtDecayBase* decayModel = EvtDecayTable::getInstance()->findDecayModel(
339 aliasInt, pythiaModeInt );
340
341 if ( decayModel != 0 ) {
342 int nModeDaug = decayModel->getNDaug();
343
344 // We need to make sure that the number of daughters match
345 if ( nDaughters == nModeDaug ) {
346 int iModeDaug( 0 );
347 for ( iModeDaug = 0; iModeDaug < nModeDaug; iModeDaug++ ) {
348 EvtId daugId = decayModel->getDaug( iModeDaug );
349 int daugPDGId = EvtPDL::getStdHep( daugId );
350 // Pythia has used the right PDG codes for this decay mode, even for conjugate modes
351 int pythiaPDGId = m_daugPDGVector[iModeDaug];
352
353 if ( daugPDGId == pythiaPDGId ) {
354 daugAliasIdVect.push_back( daugId );
355 }
356
357 } // Loop over EvtGen mode daughters
358
359 int daugAliasSize = daugAliasIdVect.size();
360 if ( daugAliasSize == nDaughters ) {
361 // All daughter Id codes are accounted for. Set the flag to stop the loop.
362 gotMode = true;
363 } else {
364 // We do not have the correct daughter ordering. Clear the id vector
365 // and try another mode.
366 daugAliasIdVect.clear();
367 }
368
369 } // Same number of daughters
370
371 } // decayModel != 0
372
373 } // Loop over available Pythia modes
374
375 if ( gotMode == false ) {
376 // We did not find a match for the daughter aliases. Just use the normal PDG codes
377 // from the Pythia decay result
378 int iPyDaug( 0 );
379 for ( iPyDaug = 0; iPyDaug < nDaughters; iPyDaug++ ) {
380 int daugPDGCode = m_daugPDGVector[iPyDaug];
381 EvtId daugPyId = EvtPDL::evtIdFromStdHep( daugPDGCode );
382 daugAliasIdVect.push_back( daugPyId );
383 }
384 }
385
386 // Make the EvtParticle daughters (with correct alias id's). Their 4-momenta are uninitialised.
387 theParent->makeDaughters( nDaughters, daugAliasIdVect );
388
389 // Now set the 4-momenta of the daughters.
390 int iDaug( 0 );
391 // Can use an iterator here, but we already had to use the vector size...
392 for ( iDaug = 0; iDaug < nDaughters; iDaug++ ) {
393 EvtParticle* theDaughter = theParent->getDaug( iDaug );
394
395 // Set the correct 4-momentum for each daughter particle.
396 if ( theDaughter != 0 ) {
397 EvtId theDaugId = daugAliasIdVect[iDaug];
398 const EvtVector4R theDaugP4 = m_daugP4Vector[iDaug];
399 theDaughter->init( theDaugId, theDaugP4 );
400 }
401 }
402}
403
404void EvtPythiaEngine::updateParticleLists()
405{
406 // Use the EvtGen decay table (decay/user.dec) to update Pythia particle
407 // decay modes. Also, make sure the generic and alias Pythia generators
408 // use the same particle data entries as defined by EvtGen (evt.pdl).
409
410 // Loop over all entries in the EvtPDL particle data table.
411 // Aliases are added at the end with id numbers equal to the
412 // original particle, but with alias integer = lastPDLEntry+1 etc..
413 int iPDL;
414 int nPDL = EvtPDL::entries();
415
416 // Reset the m_addedPDGCodes map that keeps track
417 // of any new particles added to the Pythia input data stream
418 m_addedPDGCodes.clear();
419
420 for ( iPDL = 0; iPDL < nPDL; iPDL++ ) {
421 EvtId particleId = EvtPDL::getEntry( iPDL );
422 int aliasInt = particleId.getAlias();
423
424 // Pythia and EvtGen should use the same PDG codes for particles
425 int PDGCode = EvtPDL::getStdHep( particleId );
426
427 // Update pole mass, width, lifetime and mass range
428 double mass = EvtPDL::getMeanMass( particleId );
429 double width = EvtPDL::getWidth( particleId );
430 double lifetime = EvtPDL::getctau( particleId );
431 double mmin = EvtPDL::getMinMass( particleId );
432 double mmax = EvtPDL::getMaxMass( particleId );
433
434 // Particle data pointers. The generic and alias Pythia generator pointers have
435 // their own Pythia8::ParticleData particleData public data members which contain
436 // the particle properties table. We can extract and change individual particle
437 // entries using the particleDataEntryPtr() function within ParticleData.
438 // However, we must be careful when accessing the particleData table. We must do
439 // this directly, since assigning it to another Pythia8::ParticleData object via copying
440 // or assignment will give it a different memory address and it will no longer refer to
441 // the original particleData information from the generator pointer.
442
443 ParticleDataEntryPtr entry_generic =
444 m_genericPythiaGen->particleData.particleDataEntryPtr( PDGCode );
445 ParticleDataEntryPtr entry_alias =
446 m_aliasPythiaGen->particleData.particleDataEntryPtr( PDGCode );
447
448 // Check that the PDG code is not zero/null and exclude other
449 // special cases, e.g. those reserved for internal generator use
450 if ( entry_generic != 0 && this->validPDGCode( PDGCode ) ) {
451 entry_generic->setM0( mass );
452 entry_generic->setMWidth( width );
453 entry_generic->setTau0( lifetime );
454
455 if ( std::fabs( width ) > 0.0 ) {
456 entry_generic->setMMin( mmin );
457 entry_generic->setMMax( mmax );
458 }
459 }
460
461 // Check that the PDG code is not zero/null and exclude other
462 // special cases, e.g. those reserved for internal generator use
463 if ( entry_alias != 0 && this->validPDGCode( PDGCode ) ) {
464 entry_alias->setM0( mass );
465 entry_alias->setMWidth( width );
466 entry_alias->setTau0( lifetime );
467
468 if ( std::fabs( width ) > 0.0 ) {
469 entry_alias->setMMin( mmin );
470 entry_alias->setMMax( mmax );
471 }
472 }
473
474 // Check which particles have a Pythia decay defined.
475 // Get the list of all possible decays for the particle, using the alias integer.
476 // If the particle is not actually an alias, aliasInt = idInt.
477
478 bool hasPythiaDecays = EvtDecayTable::getInstance()->hasPythia( aliasInt );
479
480 if ( hasPythiaDecays ) {
481 int isAlias = particleId.isAlias();
482
483 // Decide what generator to use depending on whether we have
484 // an aliased particle or not
485 m_thePythiaGenerator = ( isAlias == 1 ? m_aliasPythiaGen.get()
486 : m_genericPythiaGen.get() );
487
488 // Find the Pythia particle name given the standard PDG code integer
489 std::string dataName = m_thePythiaGenerator->particleData.name(
490 PDGCode );
491 bool alreadyStored = ( m_addedPDGCodes.find( abs( PDGCode ) ) !=
492 m_addedPDGCodes.end() );
493
494 if ( dataName == " " && !alreadyStored ) {
495 // Particle and its antiparticle do not exist in the Pythia database.
496 // Create a new particle, then create the new decay modes.
497 this->createPythiaParticle( particleId, PDGCode );
498 }
499
500 // For the particle, create the Pythia decay modes.
501 // Update Pythia data tables.
502 this->updatePythiaDecayTable( particleId, aliasInt, PDGCode );
503
504 } // Loop over Pythia decays
505
506 } // Loop over EvtPDL entries
507
508 //EvtGenReport(EVTGEN_INFO,"EvtGen")<<"Writing out changed generic Pythia decay list"<<endl;
509 //m_genericPythiaGen->particleData.listChanged();
510
511 //EvtGenReport(EVTGEN_INFO,"EvtGen")<<"Writing out changed alias Pythia decay list"<<endl;
512 //m_aliasPythiaGen->particleData.listChanged();
513}
514
515bool EvtPythiaEngine::validPDGCode( int PDGCode )
516{
517 // Exclude certain PDG codes: void = 0 and special values = 81 to 100, which are reserved
518 // for internal generator use (pseudoparticles) according to PDG guidelines. Also exclude
519 // nu'_tau (nu_L) = 18, which has different masses: Pythia8 = 400 GeV, EvtGen = 0 GeV.
520
521 bool isValid( true );
522
523 int absPDGCode = abs( PDGCode );
524
525 if ( absPDGCode == 0 || absPDGCode == 18 ) {
526 // Void and nu_L or nu'_tau
527 isValid = false;
528 } else if ( absPDGCode >= 81 && absPDGCode <= 100 ) {
529 // Pseudoparticles
530 isValid = false;
531 }
532
533 return isValid;
534}
535
536void EvtPythiaEngine::updatePythiaDecayTable( EvtId& particleId, int aliasInt,
537 int PDGCode )
538{
539 // Update the particle data table in Pythia.
540 // The tables store information about the allowed decay modes
541 // where the PDGId for all particles must be positive; anti-particles are stored
542 // with the corresponding particle entry.
543 // Since we do not want to implement CP violation here, just use the same branching
544 // fractions for particle and anti-particle modes.
545
546 int nModes = EvtDecayTable::getInstance()->getNModes( aliasInt );
547 int iMode( 0 );
548
549 bool firstMode( true );
550
551 // Only process positive PDG codes.
552 if ( PDGCode < 0 ) {
553 return;
554 }
555
556 // Keep track of which decay modes are Pythia decays for each aliasInt
557 std::vector<int> pythiaModes( 0 );
558
559 // Loop over the decay modes for this particle
560 for ( iMode = 0; iMode < nModes; iMode++ ) {
561 EvtDecayBase* decayModel =
562 EvtDecayTable::getInstance()->findDecayModel( aliasInt, iMode );
563
564 if ( decayModel != 0 ) {
565 int nDaug = decayModel->getNDaug();
566
567 // If the decay mode has no daughters, then that means that there will be
568 // no entries for any submode re-definitions for Pythia.
569 // This sometimes occurs for any mode using non-standard Pythia 6 codes.
570 // Do not refine the decay mode, i.e. accept the Pythia 8 default (if it exists).
571 if ( nDaug > 0 ) {
572 // Check to see if we have a Pythia decay mode
573 std::string modelName = decayModel->getModelName();
574
575 if ( modelName == "PYTHIA" ) {
576 // Keep track which decay mode is a Pythia one. We need this in order to
577 // reassign alias Id values for particles generated in the decay.
578 pythiaModes.push_back( iMode );
579
580 std::ostringstream oss;
581 oss.setf( std::ios::scientific );
582 // Write out the absolute value of the PDG code, since
583 // particles and anti-particles occupy the same part of the Pythia table.
584 oss << PDGCode;
585
586 if ( firstMode ) {
587 // Create a new channel
588 oss << ":oneChannel = ";
589 firstMode = false;
590 } else {
591 // Add the channel
592 oss << ":addChannel = ";
593 }
594
595 // Select all channels (particle and anti-particle).
596 // For CP violation, or different BFs for particle and anti-particle,
597 // use options 2 or 3 (not here).
598 int onMode( 1 );
599 oss << onMode << " ";
600
601 double BF = decayModel->getBranchingFraction();
602 oss << BF << " ";
603
604 // Need to convert the old Pythia physics mode integers with the new ones
605 // To do this, get the model argument and write a conversion method.
606 int modeInt = this->getModeInt( decayModel );
607 oss << modeInt;
608
609 int iDaug( 0 );
610 for ( iDaug = 0; iDaug < nDaug; iDaug++ ) {
611 EvtId daugId = decayModel->getDaug( iDaug );
612 int daugPDG = EvtPDL::getStdHep( daugId );
613 oss << " " << daugPDG;
614
615 } // Daughter list
616
617 m_thePythiaGenerator->readString( oss.str() );
618
619 } // is Pythia
620
621 } else {
622 EvtGenReport( EVTGEN_INFO, "EvtGen" )
623 << "Warning in EvtPythiaEngine. Trying to redefine Pythia table for "
624 << EvtPDL::name( particleId )
625 << " for a decay channel that has no daughters."
626 << " Keeping Pythia default (if available)." << endl;
627 }
628
629 } else {
630 EvtGenReport( EVTGEN_INFO, "EvtGen" )
631 << "Error in EvtPythiaEngine. decayModel is null for particle "
632 << EvtPDL::name( particleId ) << " mode number " << iMode
633 << endl;
634 }
635
636 } // Loop over modes
637
638 m_pythiaModeMap[aliasInt] = std::move(pythiaModes);
639
640 // Now, renormalise the decay branching fractions to sum to 1.0
641 std::ostringstream rescaleStr;
642 rescaleStr.setf( std::ios::scientific );
643 rescaleStr << PDGCode << ":rescaleBR = 1.0";
644
645 m_thePythiaGenerator->readString( rescaleStr.str() );
646}
647
648int EvtPythiaEngine::getModeInt( EvtDecayBase* decayModel )
649{
650 int tmpModeInt( 0 ), modeInt( 0 );
651
652 if ( decayModel != 0 ) {
653 int nVars = decayModel->getNArg();
654 // Just read the first integer, which specifies the Pythia decay model.
655 // Ignore any other values.
656 if ( nVars > 0 ) {
657 tmpModeInt = static_cast<int>( decayModel->getArg( 0 ) );
658 }
659 }
660
661 if ( m_convertPhysCodes ) {
662 // Extra code to convert the old Pythia decay model integer MDME(ICC,2) to the new one.
663 // This should be removed eventually after updating decay.dec files to use
664 // the new convention.
665
666 if ( tmpModeInt == 0 ) {
667 modeInt = 0; // phase-space
668 } else if ( tmpModeInt == 1 ) {
669 modeInt = 1; // omega or phi -> 3pi
670 } else if ( tmpModeInt == 2 ) {
671 modeInt = 11; // Dalitz decay
672 } else if ( tmpModeInt == 3 ) {
673 modeInt = 2; // V -> PS PS
674 } else if ( tmpModeInt == 4 ) {
675 modeInt = 92; // onium -> ggg or gg gamma
676 } else if ( tmpModeInt == 11 ) {
677 modeInt = 42; // phase-space of hadrons from available quarks
678 } else if ( tmpModeInt == 12 ) {
679 modeInt = 42; // phase-space for onia resonances
680 } else if ( tmpModeInt == 13 ) {
681 modeInt = 43; // phase-space of at least 3 hadrons
682 } else if ( tmpModeInt == 14 ) {
683 modeInt = 44; // phase-space of at least 4 hadrons
684 } else if ( tmpModeInt == 15 ) {
685 modeInt = 45; // phase-space of at least 5 hadrons
686 } else if ( tmpModeInt >= 22 && tmpModeInt <= 30 ) {
687 modeInt = tmpModeInt +
688 40; // phase space of hadrons with fixed multiplicity (modeInt - 60)
689 } else if ( tmpModeInt == 31 ) {
690 modeInt = 42; // two or more quarks phase-space; one spectactor quark
691 } else if ( tmpModeInt == 32 ) {
692 modeInt = 91; // qqbar or gg pair
693 } else if ( tmpModeInt == 33 ) {
694 modeInt = 0; // triplet q X qbar, where X = gluon or colour singlet (superfluous, since g's are created anyway)
695 } else if ( tmpModeInt == 41 ) {
696 modeInt = 21; // weak decay phase space, weighting nu_tau spectrum
697 } else if ( tmpModeInt == 42 ) {
698 modeInt = 22; // weak decay V-A matrix element
699 } else if ( tmpModeInt == 43 ) {
700 modeInt = 22; // weak decay V-A matrix element, quarks as jets (superfluous)
701 } else if ( tmpModeInt == 44 ) {
702 modeInt = 22; // weak decay V-A matrix element, parton showers (superfluous)
703 } else if ( tmpModeInt == 48 ) {
704 modeInt = 23; // weak decay V-A matrix element, at least 3 decay products
705 } else if ( tmpModeInt == 50 ) {
706 modeInt = 0; // default behaviour
707 } else if ( tmpModeInt == 51 ) {
708 modeInt = 0; // step threshold (channel switched off when mass daughters > mother mass
709 } else if ( tmpModeInt == 52 || tmpModeInt == 53 ) {
710 modeInt = 0; // beta-factor threshold
711 } else if ( tmpModeInt == 84 ) {
712 modeInt = 42; // unknown physics process - just use phase-space
713 } else if ( tmpModeInt == 101 ) {
714 modeInt = 0; // continuation line
715 } else if ( tmpModeInt == 102 ) {
716 modeInt = 0; // off mass shell particles.
717 } else {
718 EvtGenReport( EVTGEN_INFO, "EvtGen" )
719 << "Pythia mode integer " << tmpModeInt
720 << " is not recognised. Using phase-space model" << endl;
721 modeInt = 0; // Use phase-space for anything else
722 }
723
724 } else {
725 // No need to convert the physics mode integer code
726 modeInt = tmpModeInt;
727 }
728
729 return modeInt;
730}
731
732void EvtPythiaEngine::createPythiaParticle( EvtId& particleId, int PDGCode )
733{
734 // Use the EvtGen name, PDGId and other variables to define the new Pythia particle.
735 EvtId antiPartId = EvtPDL::chargeConj( particleId );
736
737 std::string aliasName = EvtPDL::name(
738 particleId ); // If not an alias, aliasName = normal name
739 std::string antiName = EvtPDL::name( antiPartId );
740
741 EvtSpinType::spintype spinType = EvtPDL::getSpinType( particleId );
742 int spin = EvtSpinType::getSpin2( spinType );
743
744 int charge = EvtPDL::chg3( particleId );
745
746 // Must set the correct colour type manually here, since the evt.pdl file
747 // does not store this information. This is required for quarks otherwise
748 // Pythia cannot generate the decay properly.
749 int PDGId = EvtPDL::getStdHep( particleId );
750 int colour( 0 );
751 if ( PDGId == 21 ) {
752 colour = 2; // gluons
753 } else if ( PDGId <= 8 && PDGId > 0 ) {
754 colour = 1; // single quarks
755 }
756
757 double m0 = EvtPDL::getMeanMass( particleId );
758 double mWidth = EvtPDL::getWidth( particleId );
759 double mMin = EvtPDL::getMinMass( particleId );
760 double mMax = EvtPDL::getMaxMass( particleId );
761
762 double tau0 = EvtPDL::getctau( particleId );
763
764 std::ostringstream oss;
765 oss.setf( std::ios::scientific );
766 int absPDGCode = abs( PDGCode );
767 oss << absPDGCode << ":new = " << aliasName << " " << antiName << " "
768 << spin << " " << charge << " " << colour << " " << m0 << " " << mWidth
769 << " " << mMin << " " << mMax << " " << tau0;
770
771 // Pass this information to Pythia
772 m_thePythiaGenerator->readString( oss.str() );
773
774 // Also store the absolute value of the PDG entry
775 // to keep track of which new particles have been added,
776 // which also automatically includes the anti-particle.
777 // We need to avoid creating new anti-particles when
778 // they already exist when the particle was added.
779 m_addedPDGCodes[absPDGCode] = 1;
780}
781
782void EvtPythiaEngine::updatePhysicsParameters()
783{
784 // Update any more Pythia physics (or special particle) requirements/cuts etc..
785 // This should be used if any of the Pythia 6 parameters like JetSetPar MSTJ(i) = x
786 // are needed. Such commands will need to be implemented using the new interface
787 // pythiaGenerator->readString(cmd); Here cmd is a string telling Pythia 8
788 // what physics parameters to change. This will need to be done for the generic and
789 // alias generator pointers, as appropriate.
790
791 // Set the multiplicity level for hadronic weak decays
792 std::string multiWeakCut( "ParticleDecays:multIncreaseWeak = 2.0" );
793 m_genericPythiaGen->readString( multiWeakCut );
794 m_aliasPythiaGen->readString( std::move(multiWeakCut) );
795
796 // Set the multiplicity level for all other decays
797 std::string multiCut( "ParticleDecays:multIncrease = 4.5" );
798 m_genericPythiaGen->readString( multiCut );
799 m_aliasPythiaGen->readString( std::move(multiCut) );
800
801 //Now read in any custom configuration entered in the XML
802 GeneratorCommands commands =
803 EvtExtGeneratorCommandsTable::getInstance()->getCommands( "PYTHIA" );
804 GeneratorCommands::iterator it = commands.begin();
805
806 for ( ; it != commands.end(); ++it ) {
807 Command command = *it;
808 std::vector<std::string> commandStrings;
809
810 if ( command["VERSION"] == "PYTHIA6" ) {
811 EvtGenReport( EVTGEN_INFO, "EvtGen" )
812 << "Converting Pythia 6 command: " << command["MODULE"] << "("
813 << command["PARAM"] << ")=" << command["VALUE"] << "..." << endl;
814 commandStrings = convertPythia6Command( command );
815 } else if ( command["VERSION"] == "PYTHIA8" ) {
816 commandStrings.push_back( command["MODULE"] + ":" + command["PARAM"] +
817 " = " + command["VALUE"] );
818 } else {
819 EvtGenReport( EVTGEN_ERROR, "EvtGen" )
820 << "Pythia command received by EvtPythiaEngine has bad version:"
821 << endl;
822 EvtGenReport( EVTGEN_ERROR, "EvtGen" )
823 << "Received " << command["VERSION"]
824 << " but expected PYTHIA6 or PYTHIA8." << endl;
825 EvtGenReport( EVTGEN_ERROR, "EvtGen" )
826 << "The error is likely to be in EvtDecayTable.cpp" << endl;
827 EvtGenReport( EVTGEN_ERROR, "EvtGen" )
828 << "EvtGen will now abort." << endl;
829 ::abort();
830 }
831 std::string generator = command["GENERATOR"];
832 if ( generator == "GENERIC" || generator == "Generic" ||
833 generator == "generic" || generator == "BOTH" ||
834 generator == "Both" || generator == "both" ) {
835 std::vector<std::string>::iterator it2 = commandStrings.begin();
836 for ( ; it2 != commandStrings.end(); ++it2 ) {
837 EvtGenReport( EVTGEN_INFO, "EvtGen" )
838 << "Configuring generic Pythia generator: " << ( *it2 )
839 << endl;
840 m_genericPythiaGen->readString( *it2 );
841 }
842 }
843 if ( generator == "ALIAS" || generator == "Alias" ||
844 generator == "alias" || generator == "BOTH" ||
845 generator == "Both" || generator == "both" ) {
846 std::vector<std::string>::iterator it2 = commandStrings.begin();
847 for ( ; it2 != commandStrings.end(); ++it2 ) {
848 EvtGenReport( EVTGEN_INFO, "EvtGen" )
849 << "Configuring alias Pythia generator: " << ( *it2 )
850 << endl;
851 m_aliasPythiaGen->readString( *it2 );
852 }
853 }
854 }
855}
856
857
double spin(const T &p)
Definition AtlasPID.h:1147
double charge(const T &p)
Definition AtlasPID.h:997
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition AtlasPID.h:878
std::vector< std::string > convertPythia6Command(Command command)
Pythia8::ParticleDataEntry * ParticleDataEntryPtr
static Double_t tau0
const double width
generator
Configure Herwig7 These are the commands corresponding to what would go into the regular Herwig infil...
status
Definition merge.py:16