Back to home page

EIC code displayed by LXR

 
 

    


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

0001 // HIInfo.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 
0006 // This file contains the definition of the HIInfo, EventInfo and
0007 // HIUserHooks classes, as well as the HIUnits namespace.
0008 //
0009 // EventInfo: stores full nucleon-nucleon events with corresponding Info.
0010 // HIInfo: info about a Heavy Ion run and its produced events.
0011 // HIUserHooks: User hooks for HeavyIons models.
0012 
0013 #ifndef Pythia8_HIInfo_H
0014 #define Pythia8_HIInfo_H
0015 
0016 #include "Pythia8/Pythia.h"
0017 #include "Pythia8/HIBasics.h"
0018 #include "Pythia8/HISubCollisionModel.h"
0019 
0020 namespace Pythia8 {
0021 
0022 //==========================================================================
0023 
0024 // Class for collecting info about a Heavy Ion run and its produced
0025 // events.
0026 
0027 class HIInfo {
0028 
0029 public:
0030 
0031   friend class HeavyIons;
0032   friend class Angantyr;
0033 
0034   // Constructor.
0035   HIInfo()
0036     : idProjSave(0), idTargSave(0), bSave(0.0), NSave(0), NAccSave(0),
0037       sigmaTotSave(0.0), sigmaNDSave(0.0), sigErr2TotSave(0.0),
0038       sigErr2NDSave(0.0), avNDbSave(0.0), weightSave(0.0), weightSumSave(0.0),
0039       nCollSave(10, 0), nProjSave(10, 0), nTargSave(10, 0), nFailSave(0),
0040       subCollisionsPtrSave(nullptr) {}
0041 
0042   // The impact-parameter distance in the current event.
0043   double b() const {
0044     return bSave;
0045   }
0046 
0047   // The impact parameter angle.
0048   double phi() const {
0049     return phiSave;
0050   }
0051 
0052   // The Monte Carlo integrated total cross section in the current run.
0053   double sigmaTot() const {
0054     return sigmaTotSave/millibarn;
0055   }
0056 
0057   // The estimated statistical error on sigmaTot().
0058   double sigmaTotErr() const {
0059     return sqrt(sigErr2TotSave/max(1.0, double(NSave)))/millibarn;
0060   }
0061 
0062   // The Monte Carlo integrated non-diffractive cross section in the
0063   // current run.
0064   double sigmaND() const {
0065     return sigmaNDSave/millibarn;
0066   }
0067 
0068   // The estimated statistical error on sigmaND().
0069   double sigmaNDErr() const {
0070     return sqrt(sigErr2NDSave/max(1.0, double(NSave)));
0071   }
0072 
0073   // The average NN non-diffractive impact parameter to be used to
0074   // communicate to Pythia's MPI machinery.
0075   double avNDb() const {
0076     return avNDbSave;
0077   }
0078 
0079   // The number of attempted impact parameter points.
0080   long nAttempts() const {
0081     return NSave;
0082   }
0083 
0084   // The number of produced events.
0085   long nAccepted() const {
0086     return NAccSave;
0087   }
0088 
0089   // The total number of separate sub-collisions.
0090   int nCollTot() const { return nCollSave[0]; }
0091 
0092   // The number of separate non-diffractive sub collisions in the
0093   // current event.
0094   int nCollND() const { return nCollSave[1]; }
0095 
0096   // The total number of non-diffractive sub collisions in the current event.
0097   int nCollNDTot() const { return nProjSave[1] + nTargSave[1] - nCollSave[1]; }
0098 
0099   // The number of separate single diffractive projectile excitation
0100   // sub collisions in the current event.
0101   int nCollSDP() const { return nCollSave[2]; }
0102 
0103   // The number of separate single diffractive target excitation sub
0104   // collisions in the current event.
0105   int nCollSDT() const { return nCollSave[3]; }
0106 
0107   // The number of separate double diffractive sub collisions in the
0108   // current event.
0109   int nCollDD() const { return nCollSave[4]; }
0110 
0111   // The number of separate central diffractive sub collisions in the
0112   // current event.
0113   int nCollCD() const { return nCollSave[5]; }
0114 
0115   // The number of separate elastic sub collisions.
0116   int nCollEL() const { return nCollSave[6]; }
0117 
0118   // The number of interacting projectile nucleons in the current event.
0119   int nPartProj() const { return nProjSave[0]; }
0120 
0121   // The number of absorptively wounded projectile nucleons in the
0122   // current event.
0123   int nAbsProj() const { return nProjSave[1]; }
0124 
0125   // The number of diffractively wounded projectile nucleons in the
0126   // current event.
0127   int nDiffProj() const { return nProjSave[2]; }
0128 
0129   // The number of elastically scattered projectile nucleons in the
0130   // current event.
0131   int nElProj() const { return nProjSave[3]; }
0132 
0133   // The number of interacting projectile nucleons in the current
0134   // event.
0135   int nPartTarg() const { return nTargSave[0]; }
0136 
0137   // The number of absorptively wounded projectile nucleons in the
0138   // current event.
0139   int nAbsTarg() const { return nTargSave[1]; }
0140 
0141   // The number of diffractively wounded projectile nucleons in the
0142   // current event.
0143   int nDiffTarg() const { return nTargSave[2]; }
0144 
0145   // The number of elastically scattered projectile nucleons in the
0146   // current event.
0147   int nElTarg() const { return nTargSave[3]; }
0148 
0149   // The weight for this collision.
0150   double weight() const { return weightSave; }
0151 
0152   // The sum of weights of the produced events.
0153   double weightSum() const { return weightSumSave; }
0154 
0155   // The number of failed nucleon excitations in the current event.
0156   int nFail() const {
0157     return nFailSave;
0158   }
0159 
0160   // Register a failed nucleon excitation.
0161   void failedExcitation() {
0162     ++nFailSave;
0163   }
0164 
0165 private:
0166 
0167   // Register a tried impact parameter point giving the total elastic
0168   // amplitude, the impact parameter and impact parameter generation weight.
0169   void addAttempt(double T, double bin, double phiin, double bweight);
0170 
0171   // Reweight event for whatever reason.
0172   void reweight(double w) {
0173     weightSave *= w;
0174   }
0175 
0176   // Select the primary process.
0177   void select(Info & in) {
0178     primInfo = in;
0179     primInfo.hiInfo = this;
0180   }
0181 
0182   // Accept an event and update statistics in info.
0183   void accept();
0184 
0185   // Reject an attmpted event.
0186   void reject() {}
0187 
0188   // Register a full sub collision.
0189   int addSubCollision(const SubCollision & c);
0190 
0191   // Register a participating projectile/target nucleon.
0192   int addProjectileNucleon(const Nucleon & n);
0193   int addTargetNucleon(const Nucleon & n);
0194 
0195   // Id of the colliding nuclei.
0196   int idProjSave, idTargSave;
0197 
0198   // Impact parameter.
0199   double bSave;
0200   double phiSave;
0201 
0202   // Cross section estimates.
0203   long NSave, NAccSave;
0204   double sigmaTotSave, sigmaNDSave, sigErr2TotSave, sigErr2NDSave;
0205   double avNDbSave;
0206   double weightSave, weightSumSave;
0207 
0208   // Number of collisions and paricipants. See accessor functions for
0209   // indices.
0210   vector<int> nCollSave, nProjSave, nTargSave;
0211 
0212   // Map of primary processes and the number of events and the sum of
0213   // weights.
0214   map<int,double> sumPrimW, sumPrimW2;
0215   map<int,int> NPrim;
0216   map<int,string> NamePrim;
0217 
0218   // The info object of the primary process.
0219   Info primInfo;
0220 
0221   // Number of failed nucleon excitations.
0222   int nFailSave;
0223 
0224 public:
0225 
0226   // Access to subcollision to be extracted by the user.
0227   const SubCollisionSet* subCollisionsPtr() { return subCollisionsPtrSave; }
0228 
0229 private:
0230 
0231   // Set the subcollision pointer.
0232   void setSubCollisions(const SubCollisionSet* subCollisionsPtrIn) {
0233     subCollisionsPtrSave = subCollisionsPtrIn; }
0234 
0235   // Full information about the Glauber calculation, consisting of
0236   // all subcollisions.
0237   const SubCollisionSet* subCollisionsPtrSave;
0238 
0239 };
0240 
0241 //==========================================================================
0242 
0243 // This is the heavy ion user hooks class which in the future may be
0244 // used inside a Pythia object to generate heavy ion collisons. For
0245 // now it is used outside Pythia and requires access to a number of
0246 // Pythia objects.
0247 
0248 class HIUserHooks {
0249 
0250 public:
0251 
0252   // The default constructor is empty.
0253   HIUserHooks(): idProjSave(0), idTargSave(0) {}
0254 
0255   // Virtual destructor.
0256   virtual ~HIUserHooks() {}
0257 
0258   // Initialize this user hook.
0259   virtual void init(int idProjIn, int idTargIn) {
0260     idProjSave = idProjIn;
0261     idTargSave = idTargIn;
0262   }
0263 
0264   // A user-supplied impact parameter generator.
0265   virtual bool hasImpactParameterGenerator() const { return false; }
0266   virtual shared_ptr<ImpactParameterGenerator> impactParameterGenerator()
0267     const { return nullptr; }
0268 
0269   // A user-supplied NucleusModel for the projectile and target.
0270   virtual bool hasProjectileModel() const { return false; }
0271   virtual shared_ptr<NucleusModel> projectileModel() const { return nullptr; }
0272   virtual bool hasTargetModel() const { return false; }
0273   virtual shared_ptr<NucleusModel> targetModel() const { return nullptr; }
0274 
0275   // A user-supplied SubCollisionModel for generating nucleon-nucleon
0276   // subcollisions.
0277   virtual bool hasSubCollisionModel() { return false; }
0278   virtual shared_ptr<SubCollisionModel> subCollisionModel() { return nullptr; }
0279 
0280   // A user-supplied ordering of events in (inverse) hardness.
0281   virtual bool hasEventOrdering() const { return false; }
0282   virtual double eventOrdering(const Event &, const Info &) { return -1; }
0283 
0284   // A user-supplied method for fixing up proton-neutron mismatch in
0285   // generated beams.
0286   virtual bool canFixIsoSpin() const { return false; }
0287   virtual bool fixIsoSpin(EventInfo &) { return false; }
0288 
0289   // A user-supplied method for shifting the event in impact parameter space.
0290   virtual bool canShiftEvent() const { return false; }
0291   virtual EventInfo & shiftEvent(EventInfo & ei) const { return ei; }
0292 
0293   // A user-supplied method of adding a diffractive excitation event
0294   // to another event, optionally connecting their colours.
0295   bool canAddNucleonExcitation() const { return false; }
0296   bool addNucleonExcitation(EventInfo &, EventInfo &, bool) const {
0297     return false; }
0298 
0299   // A user supplied wrapper around the Pythia::forceHadronLevel()
0300   virtual bool canForceHadronLevel() const { return false; }
0301   virtual bool forceHadronLevel(Pythia &) { return false; }
0302 
0303   // A user-supplied way of finding the remnants of an
0304   // non-diffrcative pp collision (on the target side if tside is
0305   // true) to be used to give momentum when adding.
0306   virtual bool canFindRecoilers() const { return false; }
0307   virtual vector<int>
0308   findRecoilers(const Event &, bool /* tside */, int /* beam */, int /* end */,
0309                const Vec4 & /* pdiff */, const Vec4 & /* pbeam */) const {
0310     return vector<int>();
0311   }
0312 
0313 protected:
0314 
0315   // Information set in the init() function.
0316   // The PDG id of the projectile and target nuclei.
0317   int idProjSave, idTargSave;
0318 
0319 };
0320 
0321 //==========================================================================
0322 
0323 } // end namespace Pythia8
0324 
0325 #endif // Pythia8_HIInfo_H