Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-03 08:07:10

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 /// \file Hadr09.cc
0028 /// \brief Main program of the hadronic/Hadr09 example
0029 //
0030 //------------------------------------------------------------------------
0031 // This program shows how to use the class Hadronic Generator.
0032 // The class HadronicGenerator is a kind of "hadronic generator", i.e.
0033 // provides Geant4 final states (i.e. secondary particles) produced by
0034 // hadron-nuclear inelastic collisions.
0035 // Please see the class itself for more information.
0036 //
0037 // The use of the class Hadronic Generator is very simple:
0038 // the constructor needs to be invoked only once - specifying the name
0039 // of the Geant4 "physics case" to consider ("FTFP_BERT" will be
0040 // considered as default is the name is not specified) - and then one
0041 // method needs to be called at each collision, specifying the type of
0042 // collision (hadron, energy, direction, material) to be simulated.
0043 // The class HadronicGenerator is expected to work also in a
0044 // multi-threaded environment with "external" threads (i.e. threads
0045 // that are not necessarily managed by Geant4 run-manager):
0046 // each thread should have its own instance of the class.
0047 //
0048 // See the string "***LOOKHERE***" below for the setting of parameters
0049 // of this example: the "physics case", the set of possibilities from
0050 // which to sample the projectile (i.e. whether the projectile is a
0051 // hadron or an ion - in the case of hadron projectile, a list of hadrons
0052 // is possible from which to sample at each collision; in the case of
0053 // ion projectile, only one type of ion needs to be specified),
0054 // the kinetic energy of the projectile (which can be sampled within
0055 // an interval), whether the direction of the projectile is fixed or
0056 // sampled at each collision, the target material (a list of materials
0057 // is possible, from which the target material can be sampled at each
0058 // collision, and then from this target material, the target nucleus
0059 // will be chosen randomly by Geant4 itself), and whether to print out
0060 // some information or not and how frequently.
0061 // Once a well-defined type of hadron-nucleus or nucleus-nucleus
0062 // inelastic collision has been chosen, the method
0063 //   HadronicGenerator::GenerateInteraction
0064 // returns the secondaries produced by that interaction (in the form
0065 // of a G4VParticleChange object).
0066 // Some information about this final-state is printed out as an example.
0067 //
0068 // Usage:  Hadr09
0069 //------------------------------------------------------------------------
0070 
0071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0072 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0073 
0074 #include "CLHEP/Random/Randomize.h"
0075 #include "CLHEP/Random/Ranlux64Engine.h"
0076 #include "HadronicGenerator.hh"
0077 
0078 #include "G4GenericIon.hh"
0079 #include "G4HadronicParameters.hh"
0080 #include "G4IonTable.hh"
0081 #include "G4Material.hh"
0082 #include "G4NistManager.hh"
0083 #include "G4ParticleTable.hh"
0084 #include "G4PhysicalConstants.hh"
0085 #include "G4ProcessManager.hh"
0086 #include "G4SystemOfUnits.hh"
0087 #include "G4UnitsTable.hh"
0088 #include "G4VParticleChange.hh"
0089 #include "G4ios.hh"
0090 #include "globals.hh"
0091 
0092 #include <iomanip>
0093 
0094 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0095 
0096 int main(int, char**)
0097 {
0098   G4cout << "=== Test of the HadronicGenerator ===" << G4endl;
0099 
0100   // Enable light hypernuclei and anti-hypernuclei
0101   G4HadronicParameters::Instance()->SetEnableHyperNuclei(true);
0102 
0103   // See the HadronicGenerator class for the possibilities and meaning of the "physics cases".
0104   // ( In short, it is the name of the Geant4 hadronic model used for the simulation of
0105   //   the collision, with the possibility of having a transition between two models in
0106   //   a given energy interval, as in physics lists. )
0107   const G4String namePhysics = "FTFP_BERT";  //***LOOKHERE***  PHYSICS CASE
0108   // const G4String namePhysics = "FTFP_BERT_ATL";
0109   // const G4String namePhysics = "QGSP_BERT";
0110   // const G4String namePhysics = "QGSP_BIC";
0111   // const G4String namePhysics = "FTFP_INCLXX";
0112   // const G4String namePhysics = "FTFP";
0113   // const G4String namePhysics = "QGSP";
0114   // const G4String namePhysics = "BERT";
0115   // const G4String namePhysics = "BIC";
0116   // const G4String namePhysics = "IonBIC";
0117   // const G4String namePhysics = "INCL";
0118 
0119   // The kinetic energy of the projectile will be sampled randomly, with flat probability
0120   // in the interval [minEnergy, maxEnergy].
0121   G4double minEnergy = 1.0 * CLHEP::GeV;  //***LOOKHERE***  HADRON PROJECTILE MIN Ekin
0122   G4double maxEnergy = 30.0 * CLHEP::GeV;  //***LOOKHERE***  HADRON PROJECTILE MAX Ekin
0123 
0124   const G4int numCollisions = 1000;  //***LOOKHERE***  NUMBER OF COLLISIONS
0125 
0126   // Enable or disable the print out of this program: if enabled, the number of secondaries
0127   // produced in each collisions is printed out; moreover, once every "printingGap"
0128   // collisions, the list of secondaries is printed out.
0129   const G4bool isPrintingEnabled = true;  //***LOOKHERE***  PRINT OUT ON/OFF
0130   const G4int printingGap = 100;  //***LOOKHERE***  GAP IN PRINTING
0131 
0132   // Vector of Geant4 names of hadron projectiles: one of this will be sampled randomly
0133   // (with uniform probability) for each collision, when the projectile is not a generic ion
0134   // (note that the 6 light hypernuclei and anti-hypernuclei are treated here as for the
0135   // other hadrons, not as generic ions).
0136   // Note: comment out the corresponding line in order to exclude a particle.
0137   std::vector<G4String> vecProjectiles;  //***LOOKHERE***  POSSIBLE HADRON PROJECTILES
0138   vecProjectiles.push_back("pi-");
0139   // Note: vecProjectiles.push_back( "pi0" );               // Excluded because too short-lived
0140   vecProjectiles.push_back("pi+");
0141   vecProjectiles.push_back("kaon-");
0142   vecProjectiles.push_back("kaon+");
0143   vecProjectiles.push_back("kaon0L");
0144   vecProjectiles.push_back("kaon0S");
0145   // Note: vecProjectiles.push_back( "eta" );               // Excluded because too short-lived
0146   // Note: vecProjectiles.push_back( "eta_prime" );         // Excluded because too short-lived
0147   vecProjectiles.push_back("proton");
0148   vecProjectiles.push_back("neutron");
0149   vecProjectiles.push_back("deuteron");
0150   vecProjectiles.push_back("triton");
0151   vecProjectiles.push_back("He3");
0152   vecProjectiles.push_back("alpha");
0153   vecProjectiles.push_back("lambda");
0154   vecProjectiles.push_back("sigma-");
0155   // Note: vecProjectiles.push_back( "sigma0" );            // Excluded because too short-lived
0156   vecProjectiles.push_back("sigma+");
0157   vecProjectiles.push_back("xi-");
0158   vecProjectiles.push_back("xi0");
0159   vecProjectiles.push_back("omega-");
0160   vecProjectiles.push_back("anti_proton");
0161   vecProjectiles.push_back("anti_neutron");
0162   vecProjectiles.push_back("anti_lambda");
0163   vecProjectiles.push_back("anti_sigma-");
0164   // Note: vecProjectiles.push_back( "anti_sigma0" );       // Excluded because too short-lived
0165   vecProjectiles.push_back("anti_sigma+");
0166   vecProjectiles.push_back("anti_xi-");
0167   vecProjectiles.push_back("anti_xi0");
0168   vecProjectiles.push_back("anti_omega-");
0169   vecProjectiles.push_back("anti_deuteron");
0170   vecProjectiles.push_back("anti_triton");
0171   vecProjectiles.push_back("anti_He3");
0172   vecProjectiles.push_back("anti_alpha");
0173 
0174   // Only FTFP and QGSP can handle nuclear interaction of charm and bottom hadrons
0175   if (namePhysics == "FTFP_BERT" || namePhysics == "FTFP_BERT_ATL" || namePhysics == "QGSP_BERT"
0176       || namePhysics == "QGSP_BIC" || namePhysics == "FTFP" || namePhysics == "QGSP")
0177   {
0178     // Charm and bottom hadrons
0179     vecProjectiles.push_back("D+");
0180     vecProjectiles.push_back("D-");
0181     vecProjectiles.push_back("D0");
0182     vecProjectiles.push_back("anti_D0");
0183     vecProjectiles.push_back("Ds+");
0184     vecProjectiles.push_back("Ds-");
0185     // Note: vecProjectiles.push_back( "etac" );   // Excluded because too short-lived
0186     // Note: vecProjectiles.push_back( "J/psi" );  // Excluded because too short-lived
0187     vecProjectiles.push_back("B+");
0188     vecProjectiles.push_back("B-");
0189     vecProjectiles.push_back("B0");
0190     vecProjectiles.push_back("anti_B0");
0191     vecProjectiles.push_back("Bs0");
0192     vecProjectiles.push_back("anti_Bs0");
0193     vecProjectiles.push_back("Bc+");
0194     vecProjectiles.push_back("Bc-");
0195     // Note: vecProjectiles.push_back( "Upsilon" );         // Excluded because too short-lived
0196     vecProjectiles.push_back("lambda_c+");
0197     vecProjectiles.push_back("anti_lambda_c+");
0198     // Note: vecProjectiles.push_back( "sigma_c+" );        // Excluded because too short-lived
0199     // Note: vecProjectiles.push_back( "anti_sigma_c+" );   // Excluded because too short-lived
0200     // Note: vecProjectiles.push_back( "sigma_c0" );        // Excluded because too short-lived
0201     // Note: vecProjectiles.push_back( "anti_sigma_c0" );   // Excluded because too short-lived
0202     // Note: vecProjectiles.push_back( "sigma_c++" );       // Excluded because too short-lived
0203     // Note: vecProjectiles.push_back( "anti_sigma_c++" );  // Excluded because too short-lived
0204     vecProjectiles.push_back("xi_c+");
0205     vecProjectiles.push_back("anti_xi_c+");
0206     vecProjectiles.push_back("xi_c0");
0207     vecProjectiles.push_back("anti_xi_c0");
0208     vecProjectiles.push_back("omega_c0");
0209     vecProjectiles.push_back("anti_omega_c0");
0210     vecProjectiles.push_back("lambda_b");
0211     vecProjectiles.push_back("anti_lambda_b");
0212     // Note: vecProjectiles.push_back( "sigma_b+" );       // Excluded because too short-lived
0213     // Note: vecProjectiles.push_back( "anti_sigma_b+" );  // Excluded because too short-lived
0214     // Note: vecProjectiles.push_back( "sigma_b0" );       // Excluded because too short-lived
0215     // Note: vecProjectiles.push_back( "sigma_b0" );       // Excluded because too short-lived
0216     // Note: vecProjectiles.push_back( "sigma_b-" );       // Excluded because too short-lived
0217     // Note: vecProjectiles.push_back( "anti_sigma_b-" );  // Excluded because too short-lived
0218     vecProjectiles.push_back("xi_b0");
0219     vecProjectiles.push_back("anti_xi_b0");
0220     vecProjectiles.push_back("xi_b-");
0221     vecProjectiles.push_back("anti_xi_b-");
0222     vecProjectiles.push_back("omega_b-");
0223     vecProjectiles.push_back("anti_omega_b-");
0224   }
0225 
0226   // If the hadronic interactions of light hypernuclei and anti-hypernuclei
0227   // are swtiched on, then only FTFP and INCL can handle the nuclear interactions
0228   // of light hypernuclei, and only FTFP is capable of handling the nuclear
0229   // interactions of light anti-hypernuclei.
0230   if (G4HadronicParameters::Instance()->EnableHyperNuclei()) {
0231     if (namePhysics == "FTFP_BERT" || namePhysics == "FTFP_INCLXX" || namePhysics == "FTFP"
0232         || namePhysics == "INCL")
0233     {
0234       // Light hypernuclei
0235       vecProjectiles.push_back("hypertriton");
0236       vecProjectiles.push_back("hyperalpha");
0237       vecProjectiles.push_back("hyperH4");
0238       vecProjectiles.push_back("doublehyperH4");
0239       vecProjectiles.push_back("doublehyperdoubleneutron");
0240       vecProjectiles.push_back("hyperHe5");
0241     }
0242     if (namePhysics == "FTFP_BERT" || namePhysics == "FTFP_INCLXX" || namePhysics == "FTFP") {
0243       // Light anti-hypernuclei
0244       vecProjectiles.push_back("anti_hypertriton");
0245       vecProjectiles.push_back("anti_hyperalpha");
0246       vecProjectiles.push_back("anti_hyperH4");
0247       vecProjectiles.push_back("anti_doublehyperH4");
0248       vecProjectiles.push_back("anti_doublehyperdoubleneutron");
0249       vecProjectiles.push_back("anti_hyperHe5");
0250     }
0251   }
0252 
0253   G4ParticleDefinition* projectileNucleus = nullptr;
0254   G4GenericIon* gion = G4GenericIon::GenericIon();
0255   gion->SetProcessManager(new G4ProcessManager(gion));
0256   G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
0257   G4IonTable* ions = partTable->GetIonTable();
0258   partTable->SetReadiness();
0259   ions->CreateAllIon();
0260   ions->CreateAllIsomer();
0261 
0262   const G4bool isProjectileIon = false;  //***LOOKHERE***  HADRON (false) OR ION (true) PROJECTILE?
0263   if (isProjectileIon) {
0264     minEnergy = 40.0 * 13.0 * CLHEP::GeV;  //***LOOKHERE***  ION PROJECTILE MIN Ekin
0265     maxEnergy = 40.0 * 13.0 * CLHEP::GeV;  //***LOOKHERE***  ION PROJECTILE MAX Ekin
0266     G4int ionZ = 18, ionA = 40;  //***LOOKHERE***  ION PROJECTILE (Z, A)
0267     projectileNucleus = partTable->GetIonTable()->GetIon(ionZ, ionA, 0.0);
0268   }
0269 
0270   // Vector of Geant4 NIST names of materials: one of this will be sampled randomly
0271   // (with uniform probability) for each collision and used as target material.
0272   // Note: comment out the corresponding line in order to exclude a material;
0273   //       or, vice versa, add a new line to extend the list with another material.
0274   std::vector<G4String> vecMaterials;  //***LOOKHERE*** : NIST TARGET MATERIALS
0275   vecMaterials.push_back("G4_H");
0276   vecMaterials.push_back("G4_He");
0277   vecMaterials.push_back("G4_Be");
0278   vecMaterials.push_back("G4_C");
0279   vecMaterials.push_back("G4_Al");
0280   vecMaterials.push_back("G4_Si");
0281   // vecMaterials.push_back( "G4_Sc" );
0282   vecMaterials.push_back("G4_Ar");
0283   vecMaterials.push_back("G4_Fe");
0284   vecMaterials.push_back("G4_Cu");
0285   vecMaterials.push_back("G4_W");
0286   vecMaterials.push_back("G4_Pb");
0287 
0288   const G4int numProjectiles = vecProjectiles.size();
0289   const G4int numMaterials = vecMaterials.size();
0290 
0291   G4cout << G4endl << "=================  Configuration ==================" << G4endl
0292          << "Model: " << namePhysics << G4endl << "Ekin: [ " << minEnergy / CLHEP::GeV << " , "
0293          << maxEnergy / CLHEP::GeV << " ] GeV" << G4endl
0294          << "Number of collisions:  " << numCollisions << G4endl
0295          << "Number of hadron projectiles: " << numProjectiles << G4endl
0296          << "Number of materials:   " << numMaterials << G4endl
0297          << "IsIonProjectile: " << (projectileNucleus != nullptr ? "true \t" : "false")
0298          << (projectileNucleus != nullptr ? projectileNucleus->GetParticleName() : G4String(""))
0299          << G4endl << "===================================================" << G4endl << G4endl;
0300 
0301   CLHEP::Ranlux64Engine defaultEngine(1234567, 4);
0302   CLHEP::HepRandom::setTheEngine(&defaultEngine);
0303   G4int seed = time(NULL);
0304   CLHEP::HepRandom::setTheSeed(seed);
0305   G4cout << G4endl << " Initial seed = " << seed << G4endl << G4endl;
0306 
0307   // Instanciate the HadronicGenerator providing the name of the "physics case"
0308   HadronicGenerator* theHadronicGenerator = new HadronicGenerator(namePhysics);
0309   //****************************************************************************
0310 
0311   if (theHadronicGenerator == nullptr) {
0312     G4cerr << "ERROR: theHadronicGenerator is NULL !" << G4endl;
0313     return 1;
0314   }
0315   else if (!theHadronicGenerator->IsPhysicsCaseSupported()) {
0316     G4cerr << "ERROR: this physics case is NOT supported !" << G4endl;
0317     return 2;
0318   }
0319 
0320   // Loop over the collisions
0321   G4double rnd1, rnd2, rnd3, rnd4, rnd5, rnd6, normalization, projectileEnergy;
0322   G4VParticleChange* aChange = nullptr;
0323   for (G4int i = 0; i < numCollisions; ++i) {
0324     // Draw some random numbers to select the hadron-nucleus interaction:
0325     // projectile hadron, projectile kinetic energy, projectile direction, and target material.
0326     rnd1 = CLHEP::HepRandom::getTheEngine()->flat();
0327     rnd2 = CLHEP::HepRandom::getTheEngine()->flat();
0328     rnd3 = CLHEP::HepRandom::getTheEngine()->flat();
0329     rnd4 = CLHEP::HepRandom::getTheEngine()->flat();
0330     rnd5 = CLHEP::HepRandom::getTheEngine()->flat();
0331     rnd6 = CLHEP::HepRandom::getTheEngine()->flat();
0332     // Sample the projectile kinetic energy
0333     projectileEnergy = minEnergy + rnd1 * (maxEnergy - minEnergy);
0334     if (projectileEnergy <= 0.0) projectileEnergy = minEnergy;
0335     // Sample the projectile direction
0336     normalization = 1.0 / std::sqrt(rnd2 * rnd2 + rnd3 * rnd3 + rnd4 * rnd4);
0337     const G4bool isOnSmearingDirection = true;  //***LOOKHERE***
0338     G4ThreeVector aDirection = G4ThreeVector(0.0, 0.0, 1.0);  //***LOOKHERE***
0339     if (isOnSmearingDirection) {
0340       aDirection = G4ThreeVector(normalization * rnd2, normalization * rnd3, normalization * rnd4);
0341     }
0342     // Sample the projectile hadron from the vector vecProjectiles
0343     G4int index_projectile = std::trunc(rnd5 * numProjectiles);
0344     G4String nameProjectile = vecProjectiles[index_projectile];
0345     G4ParticleDefinition* projectile = partTable->FindParticle(nameProjectile);
0346     if (projectileNucleus) {
0347       nameProjectile = projectileNucleus->GetParticleName();
0348       projectile = projectileNucleus;
0349     }
0350     // Sample the target material from the vector vecMaterials
0351     // (Note: the target nucleus will be sampled by Geant4)
0352     G4int index_material = std::trunc(rnd6 * numMaterials);
0353     G4String nameMaterial = vecMaterials[index_material];
0354     G4Material* material = G4NistManager::Instance()->FindOrBuildMaterial(nameMaterial);
0355     if (material == nullptr) {
0356       G4cerr << "ERROR: Material " << nameMaterial << " is not found !" << G4endl;
0357       return 3;
0358     }
0359     if (isPrintingEnabled) {
0360       G4cout << "\t Collision " << i << " ; projectile=" << nameProjectile;
0361       if (projectileNucleus) {
0362         G4cout << " ; Ekin[MeV]/nucleon="
0363                << projectileEnergy
0364                     / static_cast<G4double>(std::abs(projectileNucleus->GetBaryonNumber()));
0365       }
0366       else {
0367         G4cout << " ; Ekin[MeV]=" << projectileEnergy;
0368       }
0369       G4cout << " ; direction=" << aDirection << " ; material=" << nameMaterial;
0370     }
0371 
0372     // Call here the "hadronic generator" to get the secondaries produced by the hadronic collision
0373     aChange = theHadronicGenerator->GenerateInteraction(
0374       projectile, projectileEnergy,
0375       /* ********************************************** */ aDirection, material);
0376 
0377     G4int nsec = aChange ? aChange->GetNumberOfSecondaries() : 0;
0378     G4bool isPrintingOfSecondariesEnabled = false;
0379     if (isPrintingEnabled) {
0380       G4cout << G4endl << "\t --> #secondaries=" << nsec
0381              << " ; impactParameter[fm]=" << theHadronicGenerator->GetImpactParameter() / fermi
0382              << " ; #projectileSpectatorNucleons="
0383              << theHadronicGenerator->GetNumberOfProjectileSpectatorNucleons()
0384              << " ; #targetSpectatorNucleons="
0385              << theHadronicGenerator->GetNumberOfTargetSpectatorNucleons()
0386              << " ; #NNcollisions=" << theHadronicGenerator->GetNumberOfNNcollisions() << G4endl;
0387       if (i % printingGap == 0) {
0388         isPrintingOfSecondariesEnabled = true;
0389         G4cout << "\t \t List of produced secondaries: " << G4endl;
0390       }
0391     }
0392     // Loop over produced secondaries and eventually print out some information.
0393     for (G4int j = 0; j < nsec; ++j) {
0394       const G4DynamicParticle* sec = aChange->GetSecondary(j)->GetDynamicParticle();
0395       if (isPrintingOfSecondariesEnabled) {
0396         G4cout << "\t \t \t j=" << j << "\t" << sec->GetDefinition()->GetParticleName()
0397                << "\t p=" << sec->Get4Momentum() << " MeV" << G4endl;
0398       }
0399       delete aChange->GetSecondary(j);
0400     }
0401     if (aChange) aChange->Clear();
0402   }
0403 
0404   G4cout << G4endl << " Final random number = " << CLHEP::HepRandom::getTheEngine()->flat()
0405          << G4endl << "=== End of test ===" << G4endl;
0406 }
0407 
0408 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......