Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:06:40

0001 // PythiaCascade.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2024 Torbjorn Sjostrand.
0003 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
0004 // Please respect the MCnet Guidelines, see GUIDELINES for details.
0005 // Author: Torbjorn Sjostrand.
0006 
0007 #ifndef Pythia8_PythiaCascade_H
0008 #define Pythia8_PythiaCascade_H
0009 
0010 #include "Pythia8/Pythia.h"
0011 
0012 namespace Pythia8 {
0013 
0014 //==========================================================================
0015 
0016 // Wrapper class for Pythia usage in cosmic ray or detector cascades.
0017 // Here it is assumed that the user wants to control the full
0018 // evolution.  The complete event record is then kept somewhere
0019 // separate from PYTHIA, and specific input on production and decay
0020 // generates the next subevent.
0021 
0022 // Intended flow:
0023 
0024 // - init sets up all generation, given a maximum energy. This may be
0025 //   time-consuming, less so if MPI initialization data can be reused.
0026 
0027 // - sigmaSetuphN prepares a possible collision for a given hadron by
0028 //   calculating the hadron-nucleon cross section (if possible).
0029 //   Moderately time-consuming.
0030 
0031 // - sigmahA uses the hN cross section calculated by sigaSetuphN and
0032 //   returns hadron-ion cross section for a collision with a specified
0033 //   nucleus.  It can be called several times to cover a mix of target
0034 //   nuclei, with minimal time usage.
0035 
0036 // - nextColl performs the hadron-nucleus collision, as a sequences of
0037 //   hadron-nucleon ones. Can be quite time-consuming.
0038 
0039 // - nextDecay can be used anytime to decay a particle. Each
0040 //   individual decay is rather fast, but there may be many of them.
0041 
0042 // - stat can be used at end of run to give a summary of error
0043 //   messages.  Negligible time usage.
0044 
0045 // - references to particleData() and rndm() can be used in the main
0046 // - program.
0047 
0048 // Note: when a hadron interacts with a medium particle, the latter is
0049 // added to the event record. This program uses these additional
0050 // status codes for target particles in the medium:
0051 //  181: the first (or only) target nucleon in a collision.
0052 //  182: nucleons from further subcollisions with the same target nucleus.
0053 
0054 // Warning: the large boosts involved at the higher cosmic ray
0055 // energies are not perfectly managed, which leads to non-negligible
0056 // energy-momentum non-conservation. The worst (sub)events are caught
0057 // and regenerated.
0058 
0059 //--------------------------------------------------------------------------
0060 
0061 class PythiaCascade {
0062 
0063 public:
0064 
0065   // Default constructor, all setup is done in init().
0066   PythiaCascade() = default;
0067 
0068   //--------------------------------------------------------------------------
0069 
0070   // Initialize PythiaCascade for a given maximal incoming energy.
0071 
0072   // Set eMax to the maximal incoming energy (in GeV) on fixed target.
0073 
0074   // Keep listFinal = false if you want to get back the full event
0075   // record, else only the "final" particles of the collision or decay
0076   // are returned.
0077 
0078   // Keep rapidDecays = false if you want to do every decay yourself,
0079   // else all particle types with a tau0 below smallTau0 will be
0080   // decayed immediately.  The tau0 units are mm/c, so default gives c
0081   // * tau0 = 1e-10 mm = 100 fm.  Note that time dilation effects are
0082   // not included, and they can be large.
0083 
0084   // Set reuseMPI = false in case you cannot reuse the
0085   // pythiaCascade.mpi file, or any other relevant stored file. The
0086   // saved initialization must have been done with the same or a
0087   // larger eMax than the one now being used.
0088 
0089   bool init(double eMaxIn = 1e9, bool listFinalIn = false,
0090     bool rapidDecaysIn = false, double smallTau0In = 1e-10,
0091     bool reuseMPI = true, string initFile = "pythiaCascade.mpi") {
0092 
0093     // Store input for future usage.
0094     eMax        = eMaxIn;
0095     listFinal   = listFinalIn;
0096     rapidDecays = rapidDecaysIn;
0097     smallTau0   = smallTau0In;
0098 
0099     // Proton mass.
0100     mp      = pythiaMain.particleData.m0(2212);
0101 
0102     // Main Pythia object for managing the cascade evolution in a
0103     // nucleus. Can also do decays, but no hard processes.
0104     pythiaMain.readString("ProcessLevel:all = off");
0105     pythiaMain.readString("13:mayDecay  = on");
0106     pythiaMain.readString("211:mayDecay = on");
0107     pythiaMain.readString("321:mayDecay = on");
0108     pythiaMain.readString("130:mayDecay = on");
0109     pythiaMain.settings.flag("ParticleDecays:limitTau0", rapidDecays);
0110     pythiaMain.settings.parm("ParticleDecays:tau0Max", smallTau0);
0111 
0112     // Redure statistics printout to relevant ones.
0113     pythiaMain.readString("Stat:showProcessLevel = off");
0114     pythiaMain.readString("Stat:showPartonLevel = off");
0115 
0116     // Initialize. Return if failure.
0117     if (!pythiaMain.init()) return false;
0118 
0119     // Secondary Pythia object for performing individual collisions,
0120     // or decays. Variable incoming beam type and energy.
0121     pythiaColl.readString("Beams:allowVariableEnergy = on");
0122     pythiaColl.readString("Beams:allowIDAswitch = on");
0123 
0124     // Set up for fixed-target collisions.
0125     pythiaColl.readString("Beams:frameType = 3");
0126     pythiaColl.settings.parm("Beams:pzA", eMax);
0127     pythiaColl.readString("Beams:pzB = 0.");
0128 
0129     // Must use the soft and low-energy QCD processes.
0130     pythiaColl.readString("SoftQCD:all = on");
0131     pythiaColl.readString("LowEnergyQCD:all = on");
0132 
0133     // Primary (single) decay to be done by pythiaColl, to circumvent
0134     // limitTau0.
0135     pythiaColl.readString("13:mayDecay  = on");
0136     pythiaColl.readString("211:mayDecay = on");
0137     pythiaColl.readString("321:mayDecay = on");
0138     pythiaColl.readString("130:mayDecay = on");
0139 
0140     // Secondary decays to be done by pythiaMain, respecting limitTau0.
0141     pythiaColl.readString("HadronLevel:Decay = off");
0142 
0143     // Reduce printout and relax energy-momentum conservation.
0144     // (Unusually large errors unfortunate consequence of large boosts.)
0145     pythiaColl.readString("Print:quiet = on");
0146     pythiaColl.readString("Check:epTolErr = 0.01");
0147     pythiaColl.readString("Check:epTolWarn = 0.0001");
0148     pythiaColl.readString("Check:mTolErr = 0.01");
0149 
0150     // Redure statistics printout to relevant ones.
0151     pythiaColl.readString("Stat:showProcessLevel = off");
0152     pythiaColl.readString("Stat:showPartonLevel = off");
0153 
0154     // Reuse MPI initialization file if it exists; else create a new one.
0155     if (reuseMPI)
0156          pythiaColl.readString("MultipartonInteractions:reuseInit = 3");
0157     else pythiaColl.readString("MultipartonInteractions:reuseInit = 1");
0158     pythiaColl.settings.word("MultipartonInteractions:initFile", initFile);
0159 
0160     // Initialize.
0161     if (!pythiaColl.init()) return false;
0162     return true;
0163 
0164   }
0165 
0166   //--------------------------------------------------------------------------
0167 
0168   // Calculate the average number of collisions. Average number of
0169   // hadron-nucleon collisions in a hadron-nucleus one. If not in
0170   // table then interpolate, knowing that <n> - 1 propto A^{2/3}.
0171 
0172   double nCollAvg(int A) {
0173     for (int i = 0; i < nA; ++i) {
0174       if (A == tabA[i]) {
0175         return (sigmaNow < tabBorder) ? 1. + tabSlopeLo[i] * sigmaNow
0176           : 1. + tabOffset[i] + tabSlope[i] * sigmaNow;
0177       } else if (A < tabA[i]) {
0178         double nColl1 = (sigmaNow < tabBorder) ? tabSlopeLo[i - 1] * sigmaNow
0179           : tabOffset[i - 1] + tabSlope[i - 1] * sigmaNow;
0180         double nColl2 = (sigmaNow < tabBorder) ? tabSlopeLo[i] * sigmaNow
0181           : tabOffset[i] + tabSlope[i] * sigmaNow;
0182         double wt1 = (tabA[i] - A) / (tabA[i] - tabA[i - 1]);
0183         return 1. + wt1 * pow( A / tabA[i - 1], 2./3.) * nColl1
0184           + (1. - wt1) * pow( A / tabA[i], 2./3.) * nColl2;
0185       }
0186     }
0187 
0188     return numeric_limits<double>::quiet_NaN();
0189   }
0190 
0191   //--------------------------------------------------------------------------
0192 
0193   // Calculate the hadron-proton collision cross section.
0194   // Return false if not possible to find.
0195 
0196   bool sigmaSetuphN(int idNowIn, Vec4 pNowIn, double mNowIn) {
0197 
0198     // Cannot handle low-energy hadrons.
0199     if (pNowIn.e() - mNowIn < eKinMin) return false;
0200 
0201     // Cannot handle hadrons above maximum energy set at initialization.
0202     if (pNowIn.e() > eMax) {
0203       logger.ERROR_MSG("too high energy");
0204       return false;
0205     }
0206 
0207     // Save incoming quantities for reuse in later methods.
0208     idNow = idNowIn;
0209     pNow  = pNowIn;
0210     mNow  = mNowIn;
0211 
0212     // Calculate hadron-nucleon cross section. Check if cross section
0213     // vanishes.
0214     eCMNow = (pNow + Vec4(0, 0, 0, mp)).mCalc();
0215     sigmaNow = pythiaColl.getSigmaTotal(idNow, 2212, eCMNow, mNow, mp);
0216     if (sigmaNow <= 0.) {
0217       if (eCMNow - mNow - mp > eKinMin)
0218         logger.ERROR_MSG("vanishing cross section");
0219        return false;
0220     }
0221 
0222     // Done.
0223     return true;
0224 
0225   }
0226 
0227   //--------------------------------------------------------------------------
0228 
0229   // Calculate the hadron-nucleus cross section for a given nucleon
0230   // number A, using hN cross section from sigmaSetuphN. Interpolate
0231   // where not (yet) available.
0232 
0233   double sigmahA(int A) {
0234 
0235     // Restrict to allowed range 1 <= A <= 208.
0236     if (A < 1 || A > 208) {
0237       logger.ERROR_MSG("A is outside of valid range (1 <= A <= 208)");
0238       return 0.;
0239     }
0240 
0241     // Correction factor for number of h-nucleon collisions per
0242     // h-nucleus one.
0243     double sigmahA = A * sigmaNow / nCollAvg(A);
0244 
0245     // Done.
0246     return sigmahA;
0247 
0248   }
0249 
0250   //--------------------------------------------------------------------------
0251 
0252   // Generate a collision, and return the event record.
0253   // Input (Z, A) of nucleus, and optionally collision vertex.
0254 
0255   Event& nextColl(int Znow, int Anow, Vec4 vNow = Vec4() ) {
0256 
0257     // References to the two event records. Clear main event record.
0258     Event& eventMain = pythiaMain.event;
0259     Event& eventColl = pythiaColl.event;
0260     eventMain.clear();
0261 
0262     // Restrict to allowed range 1 <= A <= 208.
0263     if (Anow < 1 || Anow > 208) {
0264       logger.ERROR_MSG("A is outside of valid range (1 <= A <= 208)");
0265       return eventMain;
0266     }
0267 
0268     // Insert incoming particle in cleared main event record.
0269     eventMain.append(90,   -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
0270     int iHad = eventMain.append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
0271     eventMain[iHad].vProd(vNow);
0272 
0273     // Set up for collisions on a nucleus.
0274     int np      = Znow;
0275     int nn      = Anow - Znow;
0276     int sizeOld = 0;
0277     int sizeNew = 0;
0278     Vec4 dirNow = pNow / pNow.pAbs();
0279     Rndm& rndm  = pythiaMain.rndm;
0280 
0281     // Drop rate of geometric series. (Deuterium has corrected nCollAvg.)
0282     double probMore = 1. - 1. / nCollAvg(Anow);
0283 
0284     // Loop over varying number of hit nucleons in target nucleus.
0285     for (int iColl = 1; iColl <= Anow; ++iColl) {
0286       if (iColl > 1 && rndm.flat() > probMore) break;
0287 
0288       // Pick incoming projectile: trivial for first subcollision, else ...
0289       int iProj    = iHad;
0290       int procType = 0;
0291 
0292       // ... find highest-pLongitudinal particle from latest subcollision.
0293       if (iColl > 1) {
0294         iProj = 0;
0295         double pMax = 0.;
0296         for (int i = sizeOld; i < sizeNew; ++i)
0297         if ( eventMain[i].isFinal() && eventMain[i].isHadron()) {
0298           double pp = dot3(dirNow, eventMain[i].p());
0299           if (pp > pMax) {
0300             iProj = i;
0301             pMax  = pp;
0302           }
0303         }
0304 
0305         // No further subcollision if no particle with enough energy.
0306         if ( iProj == 0 || eventMain[iProj].e() - eventMain[iProj].m()
0307           < eKinMin) break;
0308 
0309         // Choose process; only SD or ND at perturbative energies.
0310         double eCMSub = (eventMain[iProj].p() + Vec4(0, 0, 0, mp)).mCalc();
0311         if (eCMSub > 10.) procType = (rndm.flat() < probSD) ? 4 : 1;
0312       }
0313 
0314       // Pick one p or n from target.
0315       int idProj = eventMain[iProj].id();
0316       bool doProton = rndm.flat() < (np / double(np + nn));
0317       if (doProton) np -= 1;
0318       else          nn -= 1;
0319       int idNuc = (doProton) ? 2212 : 2112;
0320 
0321       // Do a projectile-nucleon subcollision. Return empty event if failure.
0322       pythiaColl.setBeamIDs(idProj, idNuc);
0323       pythiaColl.setKinematics(eventMain[iProj].p(), Vec4());
0324       if (!pythiaColl.next(procType)) {
0325         eventMain.clear();
0326         return eventMain;
0327       }
0328 
0329       // Insert target nucleon. Mothers are (0,iProj) to mark who it
0330       // interacted with. Always use proton mass for simplicity.
0331       int statusNuc = (iColl == 1) ? -181 : -182;
0332       int iNuc = eventMain.append( idNuc, statusNuc, 0, iProj, 0, 0, 0, 0,
0333         0., 0., 0., mp, mp);
0334       eventMain[iNuc].vProdAdd(vNow);
0335 
0336       // Update full energy of the event with the proton mass.
0337       eventMain[0].e( eventMain[0].e() + mp);
0338       eventMain[0].m( eventMain[0].p().mCalc() );
0339 
0340       // Insert secondary produced particles (but skip intermediate partons)
0341       // into main event record and shift to correct production vertex.
0342       sizeOld = eventMain.size();
0343       for (int iSub = 3; iSub < eventColl.size(); ++iSub) {
0344         if (!eventColl[iSub].isFinal()) continue;
0345         int iNew = eventMain.append(eventColl[iSub]);
0346         eventMain[iNew].mothers(iNuc, iProj);
0347         eventMain[iNew].vProdAdd(vNow);
0348       }
0349       sizeNew = eventMain.size();
0350 
0351       // Update daughters of colliding hadrons and other history.
0352       eventMain[iProj].daughters(sizeOld, sizeNew - 1);
0353       eventMain[iNuc].daughters(sizeOld, sizeNew - 1);
0354       eventMain[iProj].statusNeg();
0355       eventMain[iProj].tau(0.);
0356 
0357     // End of loop over interactions in a nucleus.
0358     }
0359 
0360     // Optionally do decays of short-lived particles.
0361     if (rapidDecays) pythiaMain.moreDecays();
0362 
0363     // Optionally compress event record.
0364     if (listFinal) compress();
0365 
0366     // Return generated collision.
0367     return eventMain;
0368 
0369   }
0370 
0371   //--------------------------------------------------------------------------
0372 
0373   // Generate a particle decay, and return the event record.
0374   // You can allow sequential decays, if they occur rapidly enough.
0375 
0376   Event& nextDecay(int idNowIn, Vec4 pNowIn, double mNowIn,
0377     Vec4 vNow = Vec4() ) {
0378 
0379     // Save incoming quantities. (Not needed, but by analogy with
0380     // collisions.)
0381     idNow = idNowIn;
0382     pNow  = pNowIn;
0383     mNow  = mNowIn;
0384 
0385     // References to the event records. Clear them.
0386     Event& eventMain = pythiaMain.event;
0387     Event& eventColl = pythiaColl.event;
0388     eventMain.clear();
0389     eventColl.clear();
0390 
0391     // Insert incoming particle in cleared collision event record.
0392     eventColl.append(90,   -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
0393     int iHad = eventColl.append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
0394     eventColl[iHad].vProd(vNow);
0395 
0396     // Decay incoming particle. Return empty event if fail. Copy event record.
0397     if (!pythiaColl.moreDecays(iHad)) return eventMain;
0398     eventMain = eventColl;
0399 
0400     // Optionally do secondary decays of short-lived particles.
0401     if (rapidDecays) pythiaMain.moreDecays();
0402 
0403     // Optionally compress event record.
0404     if (listFinal) compress();
0405 
0406     // Return generated collision.
0407     return eventMain;
0408 
0409   }
0410 
0411   //--------------------------------------------------------------------------
0412 
0413   // Compress the event record by removing initial and intermediate
0414   // particles.  Keep line 0, since the += operator for event records
0415   // only adds from 1 on.
0416 
0417   void compress() {
0418 
0419     // Reference to the main event record. Original and new size.
0420     Event& eventMain = pythiaMain.event;
0421     int sizeOld = eventMain.size();
0422     int sizeNew = 1;
0423 
0424     // Loop through all particles and move up the final ones to the top.
0425     // Remove history information. Update event record size.
0426     for (int i = 1; i < sizeOld; ++i) if (eventMain[i].isFinal()) {
0427       eventMain[sizeNew] = eventMain[i];
0428       eventMain[sizeNew].mothers( 0, 0);
0429       eventMain[sizeNew].daughters( 0, 0);
0430       ++sizeNew;
0431     }
0432 
0433     // Shrink event record to new size.
0434     eventMain.popBack( sizeOld - sizeNew);
0435 
0436   }
0437 
0438   //--------------------------------------------------------------------------
0439 
0440   // Summary of aborts, errors and warnings.
0441 
0442   void stat() {
0443     pythiaMain.stat();
0444     pythiaColl.stat();
0445     logger.errorStatistics();
0446   }
0447 
0448   //--------------------------------------------------------------------------
0449 
0450   // Possibility to access particle data and random numbers from
0451   // pythiaMain.
0452 
0453   ParticleData& particleData() {return pythiaMain.particleData;}
0454 
0455   Rndm& rndm() {return pythiaMain.rndm;}
0456 
0457 //--------------------------------------------------------------------------
0458 
0459 private:
0460 
0461   // The Pythia instances used for decays and for collisions. Could
0462   // be made public, but cleaner to allow only limited access as
0463   // above.
0464   Pythia pythiaMain, pythiaColl;
0465 
0466   // Logger instance for errors in this class.
0467   Logger logger;
0468 
0469   // Save quantities.
0470   bool   listFinal, rapidDecays;
0471   int    idNow;
0472   double eMax, smallTau0, mp, mNow, eCMNow, sigmaNow;
0473   Vec4   pNow;
0474 
0475   // Fraction of single diffractive events beyond first collision in nucleus.
0476   static constexpr double probSD = 0.3;
0477 
0478   // Hadrons below the kinetic energy threshold cannot interact with medium.
0479   // (At least not using Pythia.) Used both in lab and in CM frame.
0480   static constexpr double eKinMin = 0.2;
0481 
0482   // Studied nuclei by A number, with offset and slope of <nColl>(sigma):
0483   // 1H, 2H, 4He, 9Be, 12C, 14N, 16O, 27Al, 40Ar, 56Fe, 63Cu, 84Kr,
0484   // 107Ag, 129Xe, 197Au, 208Pb.
0485   static constexpr int nA = 16;
0486   static constexpr double tabBorder = 20.;
0487   static const double tabA[nA];
0488   static const double tabOffset[nA];
0489   static const double tabSlope[nA];
0490   static const double tabSlopeLo[nA];
0491 
0492 };
0493 
0494 // Constant definitions.
0495 const double PythiaCascade::tabA[] = {1, 2, 4, 9, 12, 14, 16, 27, 40, 56,
0496   63, 84, 107, 129, 197, 208};
0497 const double PythiaCascade::tabOffset[] = {0., 0.03, 0.08, 0.15, 0.20, 0.20,
0498   0.20, 0.26, 0.30, 0.34, 0.40, 0.40, 0.40, 0.50, 0.50, 0.60};
0499 const double PythiaCascade::tabSlope[] = {0., 0.0016, 0.0033, 0.0075, 0.0092,
0500   0.0105, 0.012, 0.017, 0.022, 0.027, 0.028, 0.034, 0.040, 0.044, 0.055,
0501   0.055};
0502 const double PythiaCascade::tabSlopeLo[] = {0., 0.0031, 0.0073, 0.015, 0.0192,
0503   0.0205, 0.022, 0.03, 0.037, 0.044, 0.048, 0.054, 0.06, 0.069, 0.08, 0.085};
0504 
0505 //==========================================================================
0506 
0507 } // end namespace Pythia8
0508 
0509 #endif // end Pythia8_PythiaCascade_H