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