Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/Hadr09/Hadr09.cc-ION_PROJECTILE

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 //
 27 /// \file Hadr09.cc
 28 /// \brief Main program of the hadronic/Hadr09 example
 29 //
 30 //------------------------------------------------------------------------
 31 // This program shows how to use the class Hadronic Generator.
 32 // The class HadronicGenerator is a kind of "hadronic generator", i.e.
 33 // provides Geant4 final states (i.e. secondary particles) produced by
 34 // hadron-nuclear inelastic collisions.
 35 // Please see the class itself for more information.
 36 //
 37 // The use of the class Hadronic Generator is very simple:
 38 // the constructor needs to be invoked only once - specifying the name
 39 // of the Geant4 "physics case" to consider ("FTFP_BERT" will be
 40 // considered as default is the name is not specified) - and then one
 41 // method needs to be called at each collision, specifying the type of
 42 // collision (hadron, energy, direction, material) to be simulated.
 43 // The class HadronicGenerator is expected to work also in a
 44 // multi-threaded environment with "external" threads (i.e. threads
 45 // that are not necessarily managed by Geant4 run-manager):
 46 // each thread should have its own instance of the class.
 47 //
 48 // See the string "***LOOKHERE***" below for the setting of parameters
 49 // of this example: the "physics case", the set of possibilities from
 50 // which to sample the projectile (i.e. whether the projectile is a
 51 // hadron or an ion - in the case of hadron projectile, a list of hadrons
 52 // is possible from which to sample at each collision; in the case of
 53 // ion projectile, only one type of ion needs to be specified),
 54 // the kinetic energy of the projectile (which can be sampled within
 55 // an interval), whether the direction of the projectile is fixed or
 56 // sampled at each collision, the target material (a list of materials
 57 // is possible, from which the target material can be sampled at each
 58 // collision, and then from this target material, the target nucleus
 59 // will be chosen randomly by Geant4 itself), and whether to print out
 60 // some information or not and how frequently.
 61 // Once a well-defined type of hadron-nucleus or nucleus-nucleus
 62 // inelastic collision has been chosen, the method
 63 //   HadronicGenerator::GenerateInteraction
 64 // returns the secondaries produced by that interaction (in the form
 65 // of a G4VParticleChange object).
 66 // Some information about this final-state is printed out as an example.
 67 //
 68 // Usage:  Hadr09
 69 //------------------------------------------------------------------------
 70 
 71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 73 
 74 #include <iomanip>
 75 #include "globals.hh"
 76 #include "G4ios.hh"
 77 #include "G4PhysicalConstants.hh"
 78 #include "G4SystemOfUnits.hh"
 79 #include "G4Material.hh"
 80 #include "G4NistManager.hh"
 81 #include "G4VParticleChange.hh"
 82 #include "G4UnitsTable.hh"
 83 #include "G4SystemOfUnits.hh"
 84 #include "HadronicGenerator.hh"
 85 #include "G4GenericIon.hh"
 86 #include "G4ProcessManager.hh"
 87 #include "G4ParticleTable.hh"
 88 #include "G4IonTable.hh"
 89 #include "CLHEP/Random/Randomize.h" 
 90 #include "CLHEP/Random/Ranlux64Engine.h" 
 91 
 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 93 
 94 int main( int , char** ) {
 95   
 96   G4cout << "=== Test of the HadronicGenerator ===" << G4endl;
 97 
 98   // See the HadronicGenerator class for the possibilities and meaning of the "physics cases".
 99   // ( In short, it is the name of the Geant4 hadronic model used for the simulation of
100   //   the collision, with the possibility of having a transition between two models in
101   //   a given energy interval, as in physics lists. )
102   //const G4String namePhysics = "FTFP_BERT";        //***LOOKHERE***  PHYSICS CASE
103   //const G4String namePhysics = "FTFP_BERT_ATL";
104   //const G4String namePhysics = "QGSP_BERT";
105   //const G4String namePhysics = "QGSP_BIC";
106   //const G4String namePhysics = "FTFP_INCLXX";
107   const G4String namePhysics = "FTFP";
108   //const G4String namePhysics = "QGSP";
109   //const G4String namePhysics = "BERT";
110   //const G4String namePhysics = "BIC";
111   //const G4String namePhysics = "IonBIC";
112   //const G4String namePhysics = "INCL";
113   
114   // The kinetic energy of the projectile will be sampled randomly, with flat probability
115   // in the interval [minEnergy, maxEnergy].
116   G4double minEnergy =  1.0*CLHEP::GeV;  //***LOOKHERE***  HADRON PROJECTILE MIN Ekin
117   G4double maxEnergy = 30.0*CLHEP::GeV;  //***LOOKHERE***  HADRON PROJECTILE MAX Ekin
118   
119   const G4int numCollisions = 1000;  //***LOOKHERE***  NUMBER OF COLLISIONS
120 
121   // Enable or disable the print out of this program: if enabled, the number of secondaries
122   // produced in each collisions is printed out; moreover, once every "printingGap"
123   // collisions, the list of secondaries is printed out.
124   const G4bool isPrintingEnabled = true;           //***LOOKHERE***  PRINT OUT ON/OFF
125   const G4int  printingGap = 100;                  //***LOOKHERE***  GAP IN PRINTING
126   
127   // Vector of Geant4 names of hadron projectiles: one of this will be sampled randomly
128   // (with uniform probability) for each collision, when the projectile is not an ion.
129   // Note: comment out the corresponding line in order to exclude a particle.
130   std::vector< G4String > vecProjectiles;  //***LOOKHERE***  POSSIBLE HADRON PROJECTILES
131   vecProjectiles.push_back( "pi-" );
132   //Note: vecProjectiles.push_back( "pi0" );               // Excluded because too short-lived
133   vecProjectiles.push_back( "pi+" );
134   vecProjectiles.push_back( "kaon-" );
135   vecProjectiles.push_back( "kaon+" );
136   vecProjectiles.push_back( "kaon0L" );
137   vecProjectiles.push_back( "kaon0S" );
138   //Note: vecProjectiles.push_back( "eta" );               // Excluded because too short-lived
139   //Note: vecProjectiles.push_back( "eta_prime" );         // Excluded because too short-lived
140   vecProjectiles.push_back( "proton" );
141   vecProjectiles.push_back( "neutron" );
142   vecProjectiles.push_back( "deuteron" );
143   vecProjectiles.push_back( "triton" );
144   vecProjectiles.push_back( "He3" );
145   vecProjectiles.push_back( "alpha" );
146   vecProjectiles.push_back( "lambda" );
147   vecProjectiles.push_back( "sigma-" );
148   //Note: vecProjectiles.push_back( "sigma0" );            // Excluded because too short-lived
149   vecProjectiles.push_back( "sigma+" );
150   vecProjectiles.push_back( "xi-" );
151   vecProjectiles.push_back( "xi0" );
152   vecProjectiles.push_back( "omega-" );
153   vecProjectiles.push_back( "anti_proton" );
154   vecProjectiles.push_back( "anti_neutron" );
155   vecProjectiles.push_back( "anti_lambda" );
156   vecProjectiles.push_back( "anti_sigma-" );
157   //Note: vecProjectiles.push_back( "anti_sigma0" );       // Excluded because too short-lived
158   vecProjectiles.push_back( "anti_sigma+" );
159   vecProjectiles.push_back( "anti_xi-" );
160   vecProjectiles.push_back( "anti_xi0" );
161   vecProjectiles.push_back( "anti_omega-" );
162   vecProjectiles.push_back( "anti_deuteron" );
163   vecProjectiles.push_back( "anti_triton" );
164   vecProjectiles.push_back( "anti_He3" );
165   vecProjectiles.push_back( "anti_alpha" );
166   // Charm and bottom hadrons
167   vecProjectiles.push_back( "D+" );
168   vecProjectiles.push_back( "D-" );
169   vecProjectiles.push_back( "D0" );
170   vecProjectiles.push_back( "anti_D0" );
171   vecProjectiles.push_back( "Ds+" );
172   vecProjectiles.push_back( "Ds-" );
173   //Note: vecProjectiles.push_back( "etac" );              // Excluded because too short-lived
174   //Note: vecProjectiles.push_back( "J/psi" );             // Excluded because too short-lived
175   vecProjectiles.push_back( "B+" );
176   vecProjectiles.push_back( "B-" );
177   vecProjectiles.push_back( "B0" );
178   vecProjectiles.push_back( "anti_B0" );
179   vecProjectiles.push_back( "Bs0" );
180   vecProjectiles.push_back( "anti_Bs0" );
181   vecProjectiles.push_back( "Bc+" );
182   vecProjectiles.push_back( "Bc-" );
183   //Note: vecProjectiles.push_back( "Upsilon" );           // Excluded because too short-lived
184   vecProjectiles.push_back( "lambda_c+" );
185   vecProjectiles.push_back( "anti_lambda_c+" );
186   //Note: vecProjectiles.push_back( "sigma_c+" );          // Excluded because too short-lived
187   //Note: vecProjectiles.push_back( "anti_sigma_c+" );     // Excluded because too short-lived
188   //Note: vecProjectiles.push_back( "sigma_c0" );          // Excluded because too short-lived
189   //Note: vecProjectiles.push_back( "anti_sigma_c0" );     // Excluded because too short-lived
190   //Note: vecProjectiles.push_back( "sigma_c++" );         // Excluded because too short-lived
191   //Note: vecProjectiles.push_back( "anti_sigma_c++" );    // Excluded because too short-lived
192   vecProjectiles.push_back( "xi_c+" );
193   vecProjectiles.push_back( "anti_xi_c+" );
194   vecProjectiles.push_back( "xi_c0" );
195   vecProjectiles.push_back( "anti_xi_c0" );
196   vecProjectiles.push_back( "omega_c0" );
197   vecProjectiles.push_back( "anti_omega_c0" );
198   vecProjectiles.push_back( "lambda_b" );
199   vecProjectiles.push_back( "anti_lambda_b" );
200   //Note: vecProjectiles.push_back( "sigma_b+" );          // Excluded because too short-lived
201   //Note: vecProjectiles.push_back( "anti_sigma_b+" );     // Excluded because too short-lived  
202   //Note: vecProjectiles.push_back( "sigma_b0" );          // Excluded because too short-lived
203   //Note: vecProjectiles.push_back( "sigma_b0" );          // Excluded because too short-lived
204   //Note: vecProjectiles.push_back( "sigma_b-" );          // Excluded because too short-lived
205   //Note: vecProjectiles.push_back( "anti_sigma_b-" );     // Excluded because too short-lived  
206   vecProjectiles.push_back( "xi_b0" );
207   vecProjectiles.push_back( "anti_xi_b0" );
208   vecProjectiles.push_back( "xi_b-" );
209   vecProjectiles.push_back( "anti_xi_b-" );
210   vecProjectiles.push_back( "omega_b-" );
211   vecProjectiles.push_back( "anti_omega_b-" );
212 
213   G4ParticleDefinition* projectileNucleus = nullptr;
214   G4GenericIon* gion = G4GenericIon::GenericIon();
215   gion->SetProcessManager( new G4ProcessManager( gion ) );
216   G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
217   G4IonTable* ions = partTable->GetIonTable();
218   partTable->SetReadiness();
219   ions->CreateAllIon();
220   ions->CreateAllIsomer();
221   
222   const G4bool isProjectileIon = true;   //***LOOKHERE***  HADRON (false) OR ION (true) PROJECTILE ?
223   if ( isProjectileIon ) {
224     minEnergy = 40.0*13.0*CLHEP::GeV;    //***LOOKHERE***  ION PROJECTILE MIN Ekin
225     maxEnergy = 40.0*13.0*CLHEP::GeV;    //***LOOKHERE***  ION PROJECTILE MAX Ekin
226     G4int ionZ = 18, ionA = 40;          //***LOOKHERE***  ION PROJECTILE (Z, A)
227     projectileNucleus = partTable->GetIonTable()->GetIon( ionZ, ionA, 0.0 );
228   }
229 
230   // Vector of Geant4 NIST names of materials: one of this will be sampled randomly
231   // (with uniform probability) for each collision and used as target material.
232   // Note: comment out the corresponding line in order to exclude a material;
233   //       or, vice versa, add a new line to extend the list with another material.
234   std::vector< G4String > vecMaterials;  //***LOOKHERE*** : NIST TARGET MATERIALS
235   //vecMaterials.push_back( "G4_H" );
236   //vecMaterials.push_back( "G4_He" );
237   //vecMaterials.push_back( "G4_Be" );
238   //vecMaterials.push_back( "G4_C" );
239   //vecMaterials.push_back( "G4_Al" );
240   //vecMaterials.push_back( "G4_Si" );
241   vecMaterials.push_back( "G4_Sc" );
242   //vecMaterials.push_back( "G4_Ar" );
243   //vecMaterials.push_back( "G4_Fe" );
244   //vecMaterials.push_back( "G4_Cu" );
245   //vecMaterials.push_back( "G4_W" );
246   //vecMaterials.push_back( "G4_Pb" );
247   
248   const G4int numProjectiles = vecProjectiles.size();
249   const G4int numMaterials = vecMaterials.size();
250 
251   G4cout << G4endl
252          << "=================  Configuration ==================" << G4endl
253          << "Model: " << namePhysics << G4endl
254          << "Ekin: [ " << minEnergy/CLHEP::GeV << " , " << maxEnergy/CLHEP::GeV
255          << " ] GeV" << G4endl
256          << "Number of collisions:  " << numCollisions << G4endl
257          << "Number of hadron projectiles: " << numProjectiles << G4endl
258          << "Number of materials:   " << numMaterials   << G4endl
259    << "IsIonProjectile: " << ( projectileNucleus != nullptr ? "true \t" : "false" )
260          << ( projectileNucleus != nullptr ? projectileNucleus->GetParticleName() : "") << G4endl
261          << "===================================================" << G4endl
262          << G4endl;
263   
264   CLHEP::Ranlux64Engine defaultEngine( 1234567, 4 ); 
265   CLHEP::HepRandom::setTheEngine( &defaultEngine ); 
266   G4int seed = time( NULL ); 
267   CLHEP::HepRandom::setTheSeed( seed ); 
268   G4cout << G4endl << " Initial seed = " << seed << G4endl << G4endl; 
269   
270   // Instanciate the HadronicGenerator providing the name of the "physics case"
271   HadronicGenerator* theHadronicGenerator = new HadronicGenerator( namePhysics );
272   //****************************************************************************
273   
274   if ( theHadronicGenerator == nullptr ) {
275     G4cerr << "ERROR: theHadronicGenerator is NULL !" << G4endl;
276     return 1;
277   } else if ( ! theHadronicGenerator->IsPhysicsCaseSupported() ) {
278     G4cerr << "ERROR: this physics case is NOT supported !" << G4endl;
279     return 2;
280   }
281   
282   // Loop over the collisions
283   G4double rnd1, rnd2, rnd3, rnd4, rnd5, rnd6, normalization, projectileEnergy;
284   G4VParticleChange* aChange = nullptr;
285   for ( G4int i = 0; i < numCollisions; ++i ) {
286     // Draw some random numbers to select the hadron-nucleus interaction:
287     // projectile hadron, projectile kinetic energy, projectile direction, and target material.
288     rnd1 = CLHEP::HepRandom::getTheEngine()->flat(); 
289     rnd2 = CLHEP::HepRandom::getTheEngine()->flat();
290     rnd3 = CLHEP::HepRandom::getTheEngine()->flat();
291     rnd4 = CLHEP::HepRandom::getTheEngine()->flat();
292     rnd5 = CLHEP::HepRandom::getTheEngine()->flat();
293     rnd6 = CLHEP::HepRandom::getTheEngine()->flat();
294     // Sample the projectile kinetic energy
295     projectileEnergy = minEnergy + rnd1*( maxEnergy - minEnergy );
296     if ( projectileEnergy <= 0.0 ) projectileEnergy = minEnergy; 
297     // Sample the projectile direction
298     normalization = 1.0 / std::sqrt( rnd2*rnd2 + rnd3*rnd3 + rnd4*rnd4 );
299     const G4bool isOnSmearingDirection = false;                 //***LOOKHERE***  IF true THEN SMEAR DIRECTION 
300     G4ThreeVector aDirection = G4ThreeVector( 0.0, 0.0, 1.0 );  //***LOOKHERE***  ELSE USE THIS FIXED DIRECTION
301     if ( isOnSmearingDirection ) {
302       aDirection = G4ThreeVector( normalization*rnd2, normalization*rnd3, normalization*rnd4 );
303     } 
304     // Sample the projectile hadron from the vector vecProjectiles
305     G4int index_projectile = std::trunc( rnd5*numProjectiles );
306     G4String nameProjectile = vecProjectiles[ index_projectile ];
307     G4ParticleDefinition* projectile = partTable->FindParticle( nameProjectile );
308     if ( projectileNucleus ) {
309       nameProjectile = projectileNucleus->GetParticleName();
310       projectile = projectileNucleus;
311     }
312     // Sample the target material from the vector vecMaterials
313     // (Note: the target nucleus will be sampled by Geant4)
314     G4int index_material = std::trunc( rnd6*numMaterials );
315     G4String nameMaterial = vecMaterials[ index_material ];
316     G4Material* material = G4NistManager::Instance()->FindOrBuildMaterial( nameMaterial );
317     if ( material == nullptr ) {
318       G4cerr << "ERROR: Material " << nameMaterial << " is not found !" << G4endl;
319       return 3;
320     }
321     if ( isPrintingEnabled ) {
322       G4cout << "\t Collision " << i << " ; projectile=" << nameProjectile;
323       if ( projectileNucleus ) {
324         G4cout << " ; Ekin[MeV]/nucleon=" << projectileEnergy / 
325     static_cast< G4double >( std::abs( projectileNucleus->GetBaryonNumber() ) );
326       } else {
327         G4cout << " ; Ekin[MeV]=" << projectileEnergy; 
328       }
329       G4cout << " ; direction=" << aDirection << " ; material=" << nameMaterial;
330     }
331     
332     // Call here the "hadronic generator" to get the secondaries produced by the hadronic collision
333     aChange = theHadronicGenerator->GenerateInteraction( projectile, projectileEnergy,
334     /* ********************************************** */ aDirection, material );
335     
336     G4int nsec = aChange ? aChange->GetNumberOfSecondaries() : 0;
337     G4bool isPrintingOfSecondariesEnabled = false;
338     if ( isPrintingEnabled ) {
339       G4cout << G4endl << "\t --> #secondaries=" << nsec 
340              << " ; impactParameter[fm]=" << theHadronicGenerator->GetImpactParameter() / fermi
341        << " ; #projectileSpectatorNucleons=" << theHadronicGenerator->GetNumberOfProjectileSpectatorNucleons()
342        << " ; #targetSpectatorNucleons=" << theHadronicGenerator->GetNumberOfTargetSpectatorNucleons()
343        << " ; #NNcollisions=" << theHadronicGenerator->GetNumberOfNNcollisions() << G4endl;
344       if ( i % printingGap == 0 ) {
345         isPrintingOfSecondariesEnabled = true;
346         G4cout << "\t \t List of produced secondaries: " << G4endl;
347       }
348     }
349     // Loop over produced secondaries and eventually print out some information.
350     for ( G4int j = 0; j < nsec; ++j ) {
351       const G4DynamicParticle* sec = aChange->GetSecondary(j)->GetDynamicParticle();
352       if ( isPrintingOfSecondariesEnabled ) {
353         G4cout << "\t \t \t j=" << j << "\t" << sec->GetDefinition()->GetParticleName()
354                << "\t p=" << sec->Get4Momentum() << " MeV" << G4endl;
355       }
356       delete aChange->GetSecondary(j);
357     }
358     if ( aChange ) aChange->Clear();
359   }
360 
361   G4cout << G4endl << " Final random number = " << CLHEP::HepRandom::getTheEngine()->flat()
362          << G4endl << "=== End of test ===" << G4endl;
363 }
364 
365 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......